xref: /petsc/src/ts/utils/dmplexlandau/tutorials/ex1.c (revision 984ed092ecc7313e63e6aa733432dfb60d929fa3)
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 
199566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL,help));
209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL,NULL, "-dim", &dim, NULL));
21e0eea495SMark   /* Create a mesh */
229566063dSJacob Faibussowitsch   PetscCall(DMPlexLandauCreateVelocitySpace(PETSC_COMM_SELF, dim, "", &X, &J, &dm));
239566063dSJacob Faibussowitsch   PetscCall(DMSetUp(dm));
249566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X,&X_0));
259566063dSJacob Faibussowitsch   PetscCall(VecCopy(X,X_0));
269566063dSJacob Faibussowitsch   PetscCall(DMPlexLandauPrintNorms(X,0));
279566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, 0, 0.0));
289566063dSJacob Faibussowitsch   PetscCall(DMViewFromOptions(dm,NULL,"-dm_view"));
299566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(X,NULL,"-vec_view"));
30e0eea495SMark   /* Create timestepping solver context */
319566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_SELF,&ts));
329566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts,dm));
339566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts,&snes));
349566063dSJacob Faibussowitsch   PetscCall(SNESGetLineSearch(snes,&linesearch));
359566063dSJacob Faibussowitsch   PetscCall(SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC));
369566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,DMPlexLandauIFunction,NULL));
379566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,J,J,DMPlexLandauIJacobian,NULL));
389566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
399566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes,&ksp));
409566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp,&pc));
419566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
429566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,X));
439566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,X));
449566063dSJacob Faibussowitsch   PetscCall(DMPlexLandauPrintNorms(X,1));
459566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &time));
469566063dSJacob Faibussowitsch   PetscCall(DMSetOutputSequenceNumber(dm, 1, time));
479566063dSJacob Faibussowitsch   PetscCall(VecViewFromOptions(X,NULL,"-vec_view"));
489566063dSJacob Faibussowitsch   PetscCall(VecAXPY(X,-1,X_0));
49e0eea495SMark   /* clean up */
509566063dSJacob Faibussowitsch   PetscCall(DMPlexLandauDestroyVelocitySpace(&dm));
519566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
539566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X_0));
549566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
55b122ec5aSJacob Faibussowitsch   return 0;
56e0eea495SMark }
57e0eea495SMark 
58e0eea495SMark /*TEST
595dac466eSMark Adams   testset:
60*984ed092SMark Adams     requires: p4est !complex double defined(PETSC_USE_DMLANDAU_2D)
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