xref: /petsc/src/dm/impls/swarm/tutorials/ex1f90.F90 (revision ce78bad369055609e946c9d2c25ea67a45873e27)
14cd81d59SMatthew G. Knepley! Example program demonstrating projection between particle and finite element spaces
24cd81d59SMatthew G. Knepley      program DMSwarmTestProjection
3*ce78bad3SBarry Smith#include <petsc/finclude/petscdmplex.h>
4*ce78bad3SBarry Smith#include <petsc/finclude/petscdmswarm.h>
5*ce78bad3SBarry Smith#include <petsc/finclude/petscksp.h>
64cd81d59SMatthew G. Knepley      use petscdmplex
74cd81d59SMatthew G. Knepley      use petscdmswarm
84cd81d59SMatthew G. Knepley      use petscksp
94cd81d59SMatthew G. Knepley      implicit none
104cd81d59SMatthew G. Knepley
114cd81d59SMatthew G. Knepley      DM ::          dm, sw
124cd81d59SMatthew G. Knepley      PetscFE ::     fe
134cd81d59SMatthew G. Knepley      KSP ::         ksp
144cd81d59SMatthew G. Knepley      Mat ::         M_p, M
154cd81d59SMatthew G. Knepley      Vec ::         f, rho, rhs
1630602db0SMatthew G. Knepley      PetscInt ::    dim, Nc = 1, degree = 1, timestep = 0
174cd81d59SMatthew G. Knepley      PetscInt ::    Np = 100, p, field = 0, zero = 0, bs
184cd81d59SMatthew G. Knepley      PetscReal ::   time = 0.0,  norm
194cd81d59SMatthew G. Knepley      PetscBool ::   removePoints = PETSC_TRUE
204cd81d59SMatthew G. Knepley      PetscDataType :: dtype
214cd81d59SMatthew G. Knepley      PetscScalar, pointer :: coords(:)
224cd81d59SMatthew G. Knepley      PetscScalar, pointer :: wq(:)
234cd81d59SMatthew G. Knepley      PetscErrorCode :: ierr
244cd81d59SMatthew G. Knepley
25d8606c27SBarry Smith      PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr))
26d8606c27SBarry Smith      PetscCallA(DMCreate(PETSC_COMM_WORLD, dm, ierr))
27d8606c27SBarry Smith      PetscCallA(DMSetType(dm, DMPLEX, ierr))
28d8606c27SBarry Smith      PetscCallA(DMSetFromOptions(dm, ierr))
29d8606c27SBarry Smith      PetscCallA(DMGetDimension(dm, dim, ierr))
30*ce78bad3SBarry Smith      PetscCallA(DMViewFromOptions(dm, PETSC_NULL_OBJECT, '-dm_view', ierr))
314cd81d59SMatthew G. Knepley
324cd81d59SMatthew G. Knepley!     Create finite element space
33d8606c27SBarry Smith      PetscCallA(PetscFECreateLagrange(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, degree, PETSC_DETERMINE, fe, ierr))
34*ce78bad3SBarry Smith      PetscCallA(DMSetField(dm, field, PETSC_NULL_DMLABEL, PetscObjectCast(fe), ierr))
35d8606c27SBarry Smith      PetscCallA(DMCreateDS(dm, ierr))
36d8606c27SBarry Smith      PetscCallA(PetscFEDestroy(fe, ierr))
374cd81d59SMatthew G. Knepley
384cd81d59SMatthew G. Knepley!     Create particle swarm
39d8606c27SBarry Smith      PetscCallA(DMCreate(PETSC_COMM_WORLD, sw, ierr))
40d8606c27SBarry Smith      PetscCallA(DMSetType(sw, DMSWARM, ierr))
41d8606c27SBarry Smith      PetscCallA(DMSetDimension(sw, dim, ierr))
42d8606c27SBarry Smith      PetscCallA(DMSwarmSetType(sw, DMSWARM_PIC, ierr))
43d8606c27SBarry Smith      PetscCallA(DMSwarmSetCellDM(sw, dm, ierr))
44d8606c27SBarry Smith      PetscCallA(DMSwarmRegisterPetscDatatypeField(sw, 'w_q', Nc, PETSC_SCALAR, ierr))
45d8606c27SBarry Smith      PetscCallA(DMSwarmFinalizeFieldRegister(sw, ierr))
46d8606c27SBarry Smith      PetscCallA(DMSwarmSetLocalSizes(sw, Np, zero, ierr))
47d8606c27SBarry Smith      PetscCallA(DMSetFromOptions(sw, ierr))
48d8606c27SBarry Smith      PetscCallA(DMSwarmGetField(sw, 'w_q', bs, dtype, wq, ierr))
49d8606c27SBarry Smith      PetscCallA(DMSwarmGetField(sw, 'DMSwarmPIC_coor', bs, dtype, coords, ierr))
504cd81d59SMatthew G. Knepley      do p = 1, Np
514cd81d59SMatthew G. Knepley        coords(p*2-1) = -cos(dble(p)/dble(Np+1) * PETSC_PI)
524cd81d59SMatthew G. Knepley        coords(p*2-0) =  sin(dble(p)/dble(Np+1) * PETSC_PI)
534cd81d59SMatthew G. Knepley        wq(p)         = 1.0
544cd81d59SMatthew G. Knepley      end do
55d8606c27SBarry Smith      PetscCallA(DMSwarmRestoreField(sw, 'DMSwarmPIC_coor', bs, dtype, coords, ierr))
56d8606c27SBarry Smith      PetscCallA(DMSwarmRestoreField(sw, 'w_q', bs, dtype, wq, ierr))
57d8606c27SBarry Smith      PetscCallA(DMSwarmMigrate(sw, removePoints, ierr))
58d52c2f21SMatthew G. Knepley      PetscCallA(DMSwarmVectorDefineField(sw, 'w_q', ierr))
59*ce78bad3SBarry Smith      PetscCallA(DMViewFromOptions(sw, PETSC_NULL_OBJECT, '-swarm_view', ierr))
604cd81d59SMatthew G. Knepley
614cd81d59SMatthew G. Knepley!     Project particles to field
624cd81d59SMatthew G. Knepley!       This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE
63d8606c27SBarry Smith      PetscCallA(DMCreateMassMatrix(sw, dm, M_p, ierr))
64d8606c27SBarry Smith      PetscCallA(DMCreateGlobalVector(dm, rho, ierr))
65d8606c27SBarry Smith      PetscCallA(DMSwarmCreateGlobalVectorFromField(sw, 'w_q', f, ierr))
66d8606c27SBarry Smith      PetscCallA(MatMultTranspose(M_p, f, rho, ierr))
674cd81d59SMatthew G. Knepley
684cd81d59SMatthew G. Knepley!     Visualize mesh field
69d8606c27SBarry Smith      PetscCallA(DMSetOutputSequenceNumber(dm, timestep, time, ierr))
70*ce78bad3SBarry Smith      PetscCallA(PetscObjectViewFromOptions(PetscObjectCast(rho), PETSC_NULL_OBJECT, '-rho_view', ierr))
714cd81d59SMatthew G. Knepley
724cd81d59SMatthew G. Knepley!     Project field to particles
734cd81d59SMatthew G. Knepley!       This gives f_p = M_p^+ M f
74d8606c27SBarry Smith      PetscCallA(DMCreateMassMatrix(dm, dm, M, ierr))
75d8606c27SBarry Smith      PetscCallA(DMCreateGlobalVector(dm, rhs, ierr))
764cd81d59SMatthew G. Knepley      if (.false.) then
77d8606c27SBarry Smith         PetscCallA(MatMult(M, rho, rhs, ierr)) ! this is what you would do for and FE solve
784cd81d59SMatthew G. Knepley      else
79d8606c27SBarry Smith         PetscCallA(VecCopy(rho, rhs, ierr)) ! Identity: M^1 M rho
804cd81d59SMatthew G. Knepley      end if
81d8606c27SBarry Smith      PetscCallA(KSPCreate(PETSC_COMM_WORLD, ksp, ierr))
82d8606c27SBarry Smith      PetscCallA(KSPSetOptionsPrefix(ksp, 'ftop_', ierr))
83d8606c27SBarry Smith      PetscCallA(KSPSetFromOptions(ksp, ierr))
84d8606c27SBarry Smith      PetscCallA(KSPSetOperators(ksp, M_p, M_p, ierr))
85d8606c27SBarry Smith      PetscCallA(KSPSolveTranspose(ksp, rhs, f, ierr))
86d8606c27SBarry Smith      PetscCallA(KSPDestroy(ksp, ierr))
87d8606c27SBarry Smith      PetscCallA(VecDestroy(rhs, ierr))
88d8606c27SBarry Smith      PetscCallA(MatDestroy(M_p, ierr))
89d8606c27SBarry Smith      PetscCallA(MatDestroy(M, ierr))
904cd81d59SMatthew G. Knepley
914cd81d59SMatthew G. Knepley!     Visualize particle field
92d8606c27SBarry Smith      PetscCallA(DMSetOutputSequenceNumber(sw, timestep, time, ierr))
93*ce78bad3SBarry Smith      PetscCallA(PetscObjectViewFromOptions(PetscObjectCast(f), PETSC_NULL_OBJECT, '-weights_view', ierr))
94d8606c27SBarry Smith      PetscCallA(VecNorm(f,NORM_1,norm,ierr))
954cd81d59SMatthew G. Knepley      print *, 'Total number density = ', norm
964cd81d59SMatthew G. Knepley!     Cleanup
97d8606c27SBarry Smith      PetscCallA(DMSwarmDestroyGlobalVectorFromField(sw, 'w_q', f, ierr))
98d8606c27SBarry Smith      PetscCallA(VecDestroy(rho, ierr))
99d8606c27SBarry Smith      PetscCallA(DMDestroy(sw, ierr))
100d8606c27SBarry Smith      PetscCallA(DMDestroy(dm, ierr))
1014cd81d59SMatthew G. Knepley
102d8606c27SBarry Smith      PetscCallA(PetscFinalize(ierr))
1034cd81d59SMatthew G. Knepley      end program DMSwarmTestProjection
1044cd81d59SMatthew G. Knepley
1054cd81d59SMatthew G. Knepley!/*TEST
1064cd81d59SMatthew G. Knepley!  build:
107dfd57a17SPierre Jolivet!    requires: defined(PETSC_USING_F90FREEFORM) double !complex
1084cd81d59SMatthew G. Knepley!
1094cd81d59SMatthew G. Knepley!  test:
1104cd81d59SMatthew G. Knepley!    suffix: 0
1114cd81d59SMatthew G. Knepley!    requires: double
11230602db0SMatthew G. Knepley!    args: -dm_plex_simplex 0 -dm_plex_box_faces 40,20 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -ftop_ksp_type lsqr -ftop_pc_type none -dm_view -swarm_view
1134cd81d59SMatthew G. Knepley!    filter: grep -v DM_ | grep -v atomic
1144cd81d59SMatthew G. Knepley!    filter_output: grep -v atomic
1154cd81d59SMatthew G. Knepley!
1164cd81d59SMatthew G. Knepley!TEST*/
117