xref: /petsc/src/ts/utils/dmplexlandau/tutorials/ex1.c (revision 5dac466e5409d4ab3096829bf528bb053d52aa10)
1e0eea495SMark static char help[] = "Landau collision operator driver\n\n";
2e0eea495SMark 
3e0eea495SMark #include <petscts.h>
4e0eea495SMark #include <petsclandau.h>
5e0eea495SMark 
6e0eea495SMark int main(int argc, char **argv)
7e0eea495SMark {
8e0eea495SMark   DM             dm;
9e0eea495SMark   Vec            X,X_0;
10e0eea495SMark   PetscErrorCode ierr;
11e0eea495SMark   PetscInt       dim=2;
12e0eea495SMark   TS             ts;
13e0eea495SMark   Mat            J;
14e0eea495SMark   SNES           snes;
15e0eea495SMark   KSP            ksp;
16e0eea495SMark   PC             pc;
17e0eea495SMark   SNESLineSearch linesearch;
18e0eea495SMark   PetscReal      time;
19e0eea495SMark 
20e0eea495SMark   ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr;
21e0eea495SMark   ierr = PetscOptionsGetInt(NULL,NULL, "-dim", &dim, NULL);CHKERRQ(ierr);
22e0eea495SMark   /* Create a mesh */
23cefb98e8SMark Adams   ierr = LandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &dm);CHKERRQ(ierr);
24e0eea495SMark   ierr = LandauCreateMassMatrix(dm, NULL);CHKERRQ(ierr);
25e0eea495SMark   ierr = DMSetUp(dm);CHKERRQ(ierr);
26e0eea495SMark   ierr = VecDuplicate(X,&X_0);CHKERRQ(ierr);
27e0eea495SMark   ierr = VecCopy(X,X_0);CHKERRQ(ierr);
28e0eea495SMark   ierr = LandauPrintNorms(X,0);CHKERRQ(ierr);
29e0eea495SMark   ierr = DMSetOutputSequenceNumber(dm, 0, 0.0);CHKERRQ(ierr);
30e0eea495SMark   ierr = DMViewFromOptions(dm,NULL,"-dm_view");CHKERRQ(ierr);
31e0eea495SMark   ierr = VecViewFromOptions(X,NULL,"-vec_view");CHKERRQ(ierr);
32e0eea495SMark   /* Create timestepping solver context */
33e0eea495SMark   ierr = TSCreate(PETSC_COMM_SELF,&ts);CHKERRQ(ierr);
34e0eea495SMark   ierr = TSSetDM(ts,dm);CHKERRQ(ierr);
35e0eea495SMark   ierr = TSGetSNES(ts,&snes);CHKERRQ(ierr);
36e0eea495SMark   ierr = SNESGetLineSearch(snes,&linesearch);CHKERRQ(ierr);
37e0eea495SMark   ierr = SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC);CHKERRQ(ierr);
38e0eea495SMark   ierr = TSSetIFunction(ts,NULL,LandauIFunction,NULL);CHKERRQ(ierr);
39e0eea495SMark   ierr = TSSetIJacobian(ts,J,J,LandauIJacobian,NULL);CHKERRQ(ierr);
40e0eea495SMark   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
41e0eea495SMark   ierr = SNESGetKSP(snes,&ksp);CHKERRQ(ierr);
42e0eea495SMark   ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr);
43e0eea495SMark   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
44e0eea495SMark   ierr = TSSetSolution(ts,X);CHKERRQ(ierr);
45e0eea495SMark   ierr = TSSolve(ts,X);CHKERRQ(ierr);
46e0eea495SMark   ierr = LandauPrintNorms(X,1);CHKERRQ(ierr);
47e0eea495SMark   ierr = TSGetTime(ts, &time);CHKERRQ(ierr);
48e0eea495SMark   ierr = DMSetOutputSequenceNumber(dm, 1, time);CHKERRQ(ierr);
49e0eea495SMark   ierr = VecViewFromOptions(X,NULL,"-vec_view");CHKERRQ(ierr);
50e0eea495SMark   ierr = VecAXPY(X,-1,X_0);CHKERRQ(ierr);
51e0eea495SMark   /* clean up */
52e0eea495SMark   ierr = LandauDestroyVelocitySpace(&dm);CHKERRQ(ierr);
53e0eea495SMark   ierr = TSDestroy(&ts);CHKERRQ(ierr);
54e0eea495SMark   ierr = VecDestroy(&X);CHKERRQ(ierr);
55e0eea495SMark   ierr = VecDestroy(&X_0);CHKERRQ(ierr);
56e0eea495SMark   ierr = PetscFinalize();
57e0eea495SMark   return ierr;
58e0eea495SMark }
59e0eea495SMark 
60e0eea495SMark /*TEST
61*5dac466eSMark Adams   testset:
62*5dac466eSMark Adams     requires: p4est !complex double
63*5dac466eSMark Adams     output_file: output/ex1_0.out
64*5dac466eSMark Adams     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
65e0eea495SMark     test:
66*5dac466eSMark Adams       suffix: cpu
67*5dac466eSMark Adams       args: -dm_landau_device_type cpu
68*5dac466eSMark Adams     test:
69*5dac466eSMark Adams       suffix: kokkos
70*5dac466eSMark Adams       requires: kokkos_kernels
71*5dac466eSMark Adams       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos
72*5dac466eSMark Adams     test:
73*5dac466eSMark Adams       suffix: cuda
74*5dac466eSMark Adams       requires: cuda
75*5dac466eSMark Adams       args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda -mat_cusparse_use_cpu_solve
76e0eea495SMark 
77e0eea495SMark TEST*/
78