xref: /petsc/src/ts/utils/dmplexlandau/tutorials/ex1.c (revision 5f80ce2ab25dff0f4601e710601cbbcecf323266)
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;
21*5f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-dim", &dim, NULL));
22e0eea495SMark   /* Create a mesh */
23*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &dm));
24*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(dm));
25*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&X_0));
26*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(X,X_0));
27*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexLandauPrintNorms(X,0));
28*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetOutputSequenceNumber(dm, 0, 0.0));
29*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm,NULL,"-dm_view"));
30*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(X,NULL,"-vec_view"));
31e0eea495SMark   /* Create timestepping solver context */
32*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_SELF,&ts));
33*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,dm));
34*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSNES(ts,&snes));
35*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetLineSearch(snes,&linesearch));
36*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC));
37*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,DMPlexLandauIFunction,NULL));
38*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,J,J,DMPlexLandauIJacobian,NULL));
39*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
40*5f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetKSP(snes,&ksp));
41*5f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(ksp,&pc));
42*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
43*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,X));
44*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,X));
45*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexLandauPrintNorms(X,1));
46*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts, &time));
47*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetOutputSequenceNumber(dm, 1, time));
48*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(X,NULL,"-vec_view"));
49*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(X,-1,X_0));
50e0eea495SMark   /* clean up */
51*5f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexLandauDestroyVelocitySpace(&dm));
52*5f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
53*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
54*5f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X_0));
55e0eea495SMark   ierr = PetscFinalize();
56e0eea495SMark   return ierr;
57e0eea495SMark }
58e0eea495SMark 
59e0eea495SMark /*TEST
605dac466eSMark Adams   testset:
615dac466eSMark Adams     requires: p4est !complex double
625dac466eSMark Adams     output_file: output/ex1_0.out
635dac466eSMark 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
64e0eea495SMark     test:
655dac466eSMark Adams       suffix: cpu
665dac466eSMark Adams       args: -dm_landau_device_type cpu
675dac466eSMark Adams     test:
685dac466eSMark Adams       suffix: kokkos
695dac466eSMark Adams       requires: kokkos_kernels
705dac466eSMark Adams       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos
715dac466eSMark Adams     test:
725dac466eSMark Adams       suffix: cuda
735dac466eSMark Adams       requires: cuda
745dac466eSMark Adams       args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda -mat_cusparse_use_cpu_solve
75e0eea495SMark 
76e0eea495SMark TEST*/
77