xref: /petsc/src/ts/utils/dmplexlandau/tutorials/ex1.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
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   PetscInt       dim=2;
11e0eea495SMark   TS             ts;
12e0eea495SMark   Mat            J;
13e0eea495SMark   SNES           snes;
14e0eea495SMark   KSP            ksp;
15e0eea495SMark   PC             pc;
16e0eea495SMark   SNESLineSearch linesearch;
17e0eea495SMark   PetscReal      time;
18e0eea495SMark 
19*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc, &argv, NULL,help));
205f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetInt(NULL,NULL, "-dim", &dim, NULL));
21e0eea495SMark   /* Create a mesh */
225f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &dm));
235f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetUp(dm));
245f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDuplicate(X,&X_0));
255f80ce2aSJacob Faibussowitsch   CHKERRQ(VecCopy(X,X_0));
265f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexLandauPrintNorms(X,0));
275f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetOutputSequenceNumber(dm, 0, 0.0));
285f80ce2aSJacob Faibussowitsch   CHKERRQ(DMViewFromOptions(dm,NULL,"-dm_view"));
295f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(X,NULL,"-vec_view"));
30e0eea495SMark   /* Create timestepping solver context */
315f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_SELF,&ts));
325f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetDM(ts,dm));
335f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSNES(ts,&snes));
345f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetLineSearch(snes,&linesearch));
355f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC));
365f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIFunction(ts,NULL,DMPlexLandauIFunction,NULL));
375f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetIJacobian(ts,J,J,DMPlexLandauIJacobian,NULL));
385f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
395f80ce2aSJacob Faibussowitsch   CHKERRQ(SNESGetKSP(snes,&ksp));
405f80ce2aSJacob Faibussowitsch   CHKERRQ(KSPGetPC(ksp,&pc));
415f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
425f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetSolution(ts,X));
435f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,X));
445f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexLandauPrintNorms(X,1));
455f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTime(ts, &time));
465f80ce2aSJacob Faibussowitsch   CHKERRQ(DMSetOutputSequenceNumber(dm, 1, time));
475f80ce2aSJacob Faibussowitsch   CHKERRQ(VecViewFromOptions(X,NULL,"-vec_view"));
485f80ce2aSJacob Faibussowitsch   CHKERRQ(VecAXPY(X,-1,X_0));
49e0eea495SMark   /* clean up */
505f80ce2aSJacob Faibussowitsch   CHKERRQ(DMPlexLandauDestroyVelocitySpace(&dm));
515f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
525f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X));
535f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&X_0));
54*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
55*b122ec5aSJacob Faibussowitsch   return 0;
56e0eea495SMark }
57e0eea495SMark 
58e0eea495SMark /*TEST
595dac466eSMark Adams   testset:
605dac466eSMark Adams     requires: p4est !complex double
615dac466eSMark Adams     output_file: output/ex1_0.out
625dac466eSMark 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
63e0eea495SMark     test:
645dac466eSMark Adams       suffix: cpu
655dac466eSMark Adams       args: -dm_landau_device_type cpu
665dac466eSMark Adams     test:
675dac466eSMark Adams       suffix: kokkos
685dac466eSMark Adams       requires: kokkos_kernels
695dac466eSMark Adams       args: -dm_landau_device_type kokkos -dm_mat_type aijkokkos -dm_vec_type kokkos
705dac466eSMark Adams     test:
715dac466eSMark Adams       suffix: cuda
725dac466eSMark Adams       requires: cuda
735dac466eSMark Adams       args: -dm_landau_device_type cuda -dm_mat_type aijcusparse -dm_vec_type cuda -mat_cusparse_use_cpu_solve
74e0eea495SMark 
75e0eea495SMark TEST*/
76