xref: /petsc/src/ts/utils/dmplexlandau/tutorials/ex1f90.F90 (revision 8594ddcf839049ebf2aedeadda7d1e3dee42d2b6)
1e0eea495SMark! test phase space (Maxwellian) mesh construction (serial)
2e0eea495SMark!
3e0eea495SMark!
4e0eea495SMark!
5e0eea495SMark! Contributed by Mark Adams
6e0eea495SMarkprogram DMPlexTestLandauInterface
7e0eea495SMark  use petscts
8e0eea495SMark  use petscdmplex
9e0eea495SMark#include <petsc/finclude/petscts.h>
10e0eea495SMark#include <petsc/finclude/petscdmplex.h>
11e0eea495SMark  implicit none
12*8594ddcfSMark Adams  external DMPlexLandauIFunction
13*8594ddcfSMark Adams  external DMPlexLandauIJacobian
14e0eea495SMark  DM             dm
15e0eea495SMark  PetscInt       dim
16e0eea495SMark  PetscInt       ii
17e0eea495SMark  PetscErrorCode ierr
18e0eea495SMark  TS             ts
19e0eea495SMark  Vec            X,X_0
20e0eea495SMark  Mat            J
21e0eea495SMark  SNES           snes
22e0eea495SMark  KSP            ksp
23e0eea495SMark  PC             pc
24e0eea495SMark  SNESLineSearch linesearch
25e0eea495SMark  PetscReal      mone
26e0eea495SMark  PetscScalar    scalar
27e0eea495SMark  call PetscInitialize(PETSC_NULL_CHARACTER, ierr)
28e0eea495SMark  if (ierr .ne. 0) then
29e0eea495SMark     print*,'Unable to initialize PETSc'
30e0eea495SMark     stop
31e0eea495SMark  endif
32e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
33e0eea495SMark  !  Create mesh (DM), read in parameters, create and add f_0 (X)
34e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
35e0eea495SMark  dim = 2
36*8594ddcfSMark Adams  call DMPlexLandauCreateVelocitySpace(PETSC_COMM_WORLD, dim, '', X, J, dm, ierr);CHKERRA(ierr)
37e0eea495SMark  call DMSetUp(dm,ierr);CHKERRA(ierr)
38e0eea495SMark  call VecDuplicate(X,X_0,ierr);CHKERRA(ierr)
39e0eea495SMark  call VecCopy(X,X_0,ierr)
40e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
41e0eea495SMark  !  View
42e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
43e0eea495SMark  ii = 0
44*8594ddcfSMark Adams  call DMPlexLandauPrintNorms(X,ii,ierr);CHKERRA(ierr)
45e0eea495SMark  mone = 0;
46e0eea495SMark  call DMSetOutputSequenceNumber(dm, ii, mone, ierr);CHKERRA(ierr);
47e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
48e0eea495SMark  !    Create timestepping solver context
49e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
50e0eea495SMark  call TSCreate(PETSC_COMM_SELF,ts,ierr);CHKERRA(ierr)
51e8d2b73aSMark Adams  call TSSetOptionsPrefix(ts, 'ex1_', ierr);CHKERRA(ierr) ! should get this from the dm or give it to the dm
52e0eea495SMark  call TSSetDM(ts,dm,ierr);CHKERRA(ierr)
53e0eea495SMark  call TSGetSNES(ts,snes,ierr);CHKERRA(ierr)
54e8d2b73aSMark Adams  call SNESSetOptionsPrefix(snes, 'ex1_', ierr);CHKERRA(ierr) ! should get this from the dm or give it to the dm
55e0eea495SMark  call SNESGetLineSearch(snes,linesearch,ierr);CHKERRA(ierr)
56e0eea495SMark  call SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC,ierr);CHKERRA(ierr)
57*8594ddcfSMark Adams  call TSSetIFunction(ts,PETSC_NULL_VEC,DMPlexLandauIFunction,PETSC_NULL_VEC,ierr);CHKERRA(ierr)
58*8594ddcfSMark Adams  call TSSetIJacobian(ts,J,J,DMPlexLandauIJacobian,PETSC_NULL_VEC,ierr);CHKERRA(ierr)
59e0eea495SMark  call TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER,ierr);CHKERRA(ierr)
60e0eea495SMark
61e0eea495SMark  call SNESGetKSP(snes,ksp,ierr);CHKERRA(ierr)
62e8d2b73aSMark Adams  call KSPSetOptionsPrefix(ksp, 'ex1_', ierr);CHKERRA(ierr) ! should get this from the dm or give it to the dm
63e0eea495SMark  call KSPGetPC(ksp,pc,ierr);CHKERRA(ierr)
64e8d2b73aSMark Adams  call PCSetOptionsPrefix(pc, 'ex1_', ierr);CHKERRA(ierr) ! should get this from the dm or give it to the dm
65e0eea495SMark
66e0eea495SMark  call TSSetFromOptions(ts,ierr);CHKERRA(ierr)
67e0eea495SMark  call TSSetSolution(ts,X,ierr);CHKERRA(ierr)
68e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
69e0eea495SMark  !  Solve nonlinear system
70e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71e0eea495SMark  call TSSolve(ts,X,ierr);CHKERRA(ierr)
72e0eea495SMark  ii = 1
73*8594ddcfSMark Adams  call DMPlexLandauPrintNorms(X,ii,ierr);CHKERRA(ierr)
74e0eea495SMark  call TSGetTime(ts, mone, ierr);CHKERRA(ierr);
75e0eea495SMark  call DMSetOutputSequenceNumber(dm, ii, mone, ierr);CHKERRA(ierr);
76e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77e0eea495SMark  !  remove f_0
78e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
79e0eea495SMark  scalar = -1.
80e0eea495SMark  call VecAXPY(X,scalar,X_0,ierr);CHKERRA(ierr)
81*8594ddcfSMark Adams  call DMPlexLandauDestroyVelocitySpace(dm, ierr);CHKERRA(ierr)
82e0eea495SMark  call TSDestroy(ts, ierr);CHKERRA(ierr)
83e0eea495SMark  call VecDestroy(X, ierr);CHKERRA(ierr)
84e0eea495SMark  call VecDestroy(X_0, ierr);CHKERRA(ierr)
85e0eea495SMark  call PetscFinalize(ierr)
86e0eea495SMarkend program DMPlexTestLandauInterface
87e0eea495SMark
88e0eea495SMark!/*TEST
89e0eea495SMark!  build:
90dfd57a17SPierre Jolivet!    requires: defined(PETSC_USING_F90FREEFORM)
91e0eea495SMark!
92e0eea495SMark!  test:
93e0eea495SMark!    suffix: 0
948fdabdddSMark Adams!    requires: p4est !complex  !kokkos_kernels !cuda
95365b711fSMark 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 -ex1_ts_monitor -ex1_snes_rtol 1.e-14 -ex1_snes_stol 1.e-14 -ex1_snes_monitor -ex1_snes_converged_reason -ex1_ts_type arkimex -ex1_ts_arkimex_type 1bee -ex1_ts_max_snes_failures -1 -ex1_ts_rtol 1e-1 -ex1_ts_dt 1.e-1 -ex1_ts_max_time 1 -ex1_ts_adapt_clip .5,1.25 -ex1_ts_adapt_scale_solve_failed 0.75 -ex1_ts_adapt_time_step_increase_delay 5 -ex1_ts_max_steps 1 -ex1_pc_type lu -ex1_ksp_type preonly -dm_landau_amr_levels_max 2,1 -dm_landau_device_type cpu
96e0eea495SMark!
97e0eea495SMark!TEST*/
98