xref: /petsc/src/ts/utils/dmplexlandau/tutorials/ex1f90.F90 (revision 930e68a5a1dbab7612595fd12ba2e3af4e1d5d80)
1e0eea495SMark! test phase space (Maxwellian) mesh construction (serial)
2e0eea495SMark!
3e0eea495SMark!run:
4e0eea495SMark!       -${MPIEXEC} ....
5e0eea495SMark!       -@${PETSC_DIR}/lib/petsc/bin/petsc_gen_xdmf.py *.h5
6e0eea495SMark!
7e0eea495SMark!
8e0eea495SMark! Contributed by Mark Adams
9e0eea495SMarkprogram DMPlexTestLandauInterface
10e0eea495SMark  use petscts
11e0eea495SMark  use petscdmplex
12e0eea495SMark#include <petsc/finclude/petscts.h>
13e0eea495SMark#include <petsc/finclude/petscdmplex.h>
14e0eea495SMark  implicit none
15e0eea495SMark  external LandauIFunction
16e0eea495SMark  external LandauIJacobian
17e0eea495SMark  DM             dm
18e0eea495SMark  PetscInt       dim
19e0eea495SMark  PetscInt       ii
20e0eea495SMark  PetscErrorCode ierr
21e0eea495SMark  TS             ts
22e0eea495SMark  Vec            X,X_0
23e0eea495SMark  Mat            J
24e0eea495SMark  SNES           snes
25e0eea495SMark  KSP            ksp
26e0eea495SMark  PC             pc
27e0eea495SMark  SNESLineSearch linesearch
28e0eea495SMark  PetscReal      mone
29e0eea495SMark  PetscScalar    scalar
30e0eea495SMark  call PetscInitialize(PETSC_NULL_CHARACTER, ierr)
31e0eea495SMark  if (ierr .ne. 0) then
32e0eea495SMark     print*,'Unable to initialize PETSc'
33e0eea495SMark     stop
34e0eea495SMark  endif
35e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
36e0eea495SMark  !  Create mesh (DM), read in parameters, create and add f_0 (X)
37e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
38e0eea495SMark  dim = 2
39*930e68a5SMark Adams  call LandauCreateVelocitySpace(PETSC_COMM_WORLD, dim, '', X, J, dm, ierr);CHKERRA(ierr)
40e0eea495SMark  call LandauCreateMassMatrix(dm, PETSC_NULL_MAT, ierr);CHKERRA(ierr)
41e0eea495SMark  call DMSetUp(dm,ierr);CHKERRA(ierr)
42e0eea495SMark  call VecDuplicate(X,X_0,ierr);CHKERRA(ierr)
43e0eea495SMark  call VecCopy(X,X_0,ierr)
44e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
45e0eea495SMark  !  View
46e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
47e0eea495SMark  ii = 0
48e0eea495SMark  call LandauPrintNorms(X,ii,ierr);CHKERRA(ierr)
49e0eea495SMark  mone = 0;
50e0eea495SMark  call DMSetOutputSequenceNumber(dm, ii, mone, ierr);CHKERRA(ierr);
51e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
52e0eea495SMark  !    Create timestepping solver context
53e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54e0eea495SMark  call TSCreate(PETSC_COMM_SELF,ts,ierr);CHKERRA(ierr)
55e0eea495SMark  call TSSetOptionsPrefix(ts, 'ex1f90_', ierr);CHKERRA(ierr) ! should get this from the dm or give it to the dm
56e0eea495SMark  call TSSetDM(ts,dm,ierr);CHKERRA(ierr)
57e0eea495SMark  call TSGetSNES(ts,snes,ierr);CHKERRA(ierr)
58e0eea495SMark  call SNESSetOptionsPrefix(snes, 'ex1f90_', ierr);CHKERRA(ierr) ! should get this from the dm or give it to the dm
59e0eea495SMark  call SNESGetLineSearch(snes,linesearch,ierr);CHKERRA(ierr)
60e0eea495SMark  call SNESLineSearchSetType(linesearch,SNESLINESEARCHBASIC,ierr);CHKERRA(ierr)
61e0eea495SMark  call TSSetIFunction(ts,PETSC_NULL_VEC,LandauIFunction,PETSC_NULL_VEC,ierr);CHKERRA(ierr)
62e0eea495SMark  call TSSetIJacobian(ts,J,J,LandauIJacobian,PETSC_NULL_VEC,ierr);CHKERRA(ierr)
63e0eea495SMark  call TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER,ierr);CHKERRA(ierr)
64e0eea495SMark
65e0eea495SMark  call SNESGetKSP(snes,ksp,ierr);CHKERRA(ierr)
66e0eea495SMark  call KSPSetOptionsPrefix(ksp, 'ex1f90_', ierr);CHKERRA(ierr) ! should get this from the dm or give it to the dm
67e0eea495SMark  call KSPGetPC(ksp,pc,ierr);CHKERRA(ierr)
68e0eea495SMark  call PCSetOptionsPrefix(pc, 'ex1f90_', ierr);CHKERRA(ierr) ! should get this from the dm or give it to the dm
69e0eea495SMark
70e0eea495SMark  call TSSetFromOptions(ts,ierr);CHKERRA(ierr)
71e0eea495SMark  call TSSetSolution(ts,X,ierr);CHKERRA(ierr)
72e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73e0eea495SMark  !  Solve nonlinear system
74e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
75e0eea495SMark  call TSSolve(ts,X,ierr);CHKERRA(ierr)
76e0eea495SMark  ii = 1
77e0eea495SMark  call LandauPrintNorms(X,ii,ierr);CHKERRA(ierr)
78e0eea495SMark  call TSGetTime(ts, mone, ierr);CHKERRA(ierr);
79e0eea495SMark  call DMSetOutputSequenceNumber(dm, ii, mone, ierr);CHKERRA(ierr);
80e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
81e0eea495SMark  !  remove f_0
82e0eea495SMark  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
83e0eea495SMark  scalar = -1.
84e0eea495SMark  call VecAXPY(X,scalar,X_0,ierr);CHKERRA(ierr)
85e0eea495SMark  call LandauDestroyVelocitySpace(dm, ierr);CHKERRA(ierr)
86e0eea495SMark  call TSDestroy(ts, ierr);CHKERRA(ierr)
87e0eea495SMark  call VecDestroy(X, ierr);CHKERRA(ierr)
88e0eea495SMark  call VecDestroy(X_0, ierr);CHKERRA(ierr)
89e0eea495SMark  call PetscFinalize(ierr)
90e0eea495SMarkend program DMPlexTestLandauInterface
91e0eea495SMark
92e0eea495SMark!/*TEST
93e0eea495SMark!  build:
94e0eea495SMark!    requires: define(PETSC_USING_F90FREEFORM)
95e0eea495SMark!
96e0eea495SMark!  test:
97e0eea495SMark!    suffix: 0
98e0eea495SMark!    requires: p4est !complex
99a587d139SMark!    args: -petscspace_degree 3 -petscspace_poly_tensor 1 -dm_landau_type p4est -info :dm,tsadapt -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 -ex1f90_ts_monitor -ex1f90_snes_rtol 1.e-6 -ex1f90_snes_monitor -ex1f90_snes_converged_reason -ex1f90_ts_type arkimex -ex1f90_ts_arkimex_type 1bee -ex1f90_ts_max_snes_failures -1 -ex1f90_ts_rtol 1e-4 -ex1f90_ts_dt 1.e-6 -ex1f90_ts_max_time 1 -ex1f90_ts_adapt_clip .5,1.25 -ex1f90_ts_adapt_scale_solve_failed 0.75 -ex1f90_ts_adapt_time_step_increase_delay 5 -ex1f90_ts_max_steps 1 -ex1f90_pc_type lu -ex1f90_ksp_type preonly -dm_landau_amr_levels_max 7 -dm_landau_domain_radius 5 -dm_landau_amr_re_levels 0 -dm_landau_re_radius 1 -dm_landau_amr_z_refine1 1 -dm_landau_amr_z_refine2 0 -dm_landau_amr_post_refine 0 -dm_landau_z_radius1 .1 -dm_landau_z_radius2 .1 -dm_refine 1 -dm_landau_device_type cpu
100e0eea495SMark!
101e0eea495SMark!TEST*/
102