xref: /petsc/src/dm/impls/swarm/tutorials/ex1f90.F90 (revision dfd57a172ac9fa6c7b5fe6de6ab5df85cefc2996)
14cd81d59SMatthew G. Knepley! Example program demonstrating projection between particle and finite element spaces
24cd81d59SMatthew G. Knepley      program DMSwarmTestProjection
34cd81d59SMatthew G. Knepley#include "petsc/finclude/petscdmplex.h"
44cd81d59SMatthew G. Knepley#include "petsc/finclude/petscdmswarm.h"
54cd81d59SMatthew G. Knepley#include "petsc/finclude/petscksp.h"
64cd81d59SMatthew G. Knepley      use petscdmplex
74cd81d59SMatthew G. Knepley      use petscdmswarm
84cd81d59SMatthew G. Knepley      use petscdt
94cd81d59SMatthew G. Knepley      use petscksp
104cd81d59SMatthew G. Knepley      use petscsys
114cd81d59SMatthew G. Knepley      implicit none
124cd81d59SMatthew G. Knepley
134cd81d59SMatthew G. Knepley      DM ::          dm, sw
144cd81d59SMatthew G. Knepley      PetscFE ::     fe
154cd81d59SMatthew G. Knepley      KSP ::         ksp
164cd81d59SMatthew G. Knepley      Mat ::         M_p, M
174cd81d59SMatthew G. Knepley      Vec ::         f, rho, rhs
1830602db0SMatthew G. Knepley      PetscInt ::    dim, Nc = 1, degree = 1, timestep = 0
194cd81d59SMatthew G. Knepley      PetscInt ::    Np = 100, p, field = 0, zero = 0, bs
204cd81d59SMatthew G. Knepley      PetscReal ::   time = 0.0,  norm
214cd81d59SMatthew G. Knepley      PetscBool ::   removePoints = PETSC_TRUE
224cd81d59SMatthew G. Knepley      PetscDataType :: dtype
234cd81d59SMatthew G. Knepley      PetscScalar, pointer :: coords(:)
244cd81d59SMatthew G. Knepley      PetscScalar, pointer :: wq(:)
254cd81d59SMatthew G. Knepley      PetscErrorCode :: ierr
264cd81d59SMatthew G. Knepley
274cd81d59SMatthew G. Knepley      call PetscInitialize(PETSC_NULL_CHARACTER, ierr)
284cd81d59SMatthew G. Knepley      if (ierr .ne. 0) then
294cd81d59SMatthew G. Knepley        print*,'Unable to initialize PETSc'
304cd81d59SMatthew G. Knepley        stop
314cd81d59SMatthew G. Knepley      endif
324cd81d59SMatthew G. Knepley
3330602db0SMatthew G. Knepley      call DMCreate(PETSC_COMM_WORLD, dm, ierr);CHKERRA(ierr)
3430602db0SMatthew G. Knepley      call DMSetType(dm, DMPLEX, ierr);CHKERRA(ierr)
3530602db0SMatthew G. Knepley      call DMSetFromOptions(dm, ierr);CHKERRA(ierr)
3630602db0SMatthew G. Knepley      call DMGetDimension(dm, dim, ierr);CHKERRA(ierr)
374cd81d59SMatthew G. Knepley      call DMViewFromOptions(dm, PETSC_NULL_VEC, '-dm_view', ierr);CHKERRA(ierr)
384cd81d59SMatthew G. Knepley
394cd81d59SMatthew G. Knepley!     Create finite element space
404cd81d59SMatthew G. Knepley      call PetscFECreateLagrange(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, degree, PETSC_DETERMINE, fe, ierr);CHKERRA(ierr)
414cd81d59SMatthew G. Knepley      call DMSetField(dm, field, PETSC_NULL_DMLABEL, fe, ierr);CHKERRA(ierr)
424cd81d59SMatthew G. Knepley      call DMCreateDS(dm, ierr);CHKERRA(ierr)
434cd81d59SMatthew G. Knepley      call PetscFEDestroy(fe, ierr);CHKERRA(ierr)
444cd81d59SMatthew G. Knepley
454cd81d59SMatthew G. Knepley!     Create particle swarm
464cd81d59SMatthew G. Knepley      call DMCreate(PETSC_COMM_WORLD, sw, ierr);CHKERRA(ierr)
474cd81d59SMatthew G. Knepley      call DMSetType(sw, DMSWARM, ierr);CHKERRA(ierr)
484cd81d59SMatthew G. Knepley      call DMSetDimension(sw, dim, ierr);CHKERRA(ierr)
494cd81d59SMatthew G. Knepley      call DMSwarmSetType(sw, DMSWARM_PIC, ierr);CHKERRA(ierr)
504cd81d59SMatthew G. Knepley      call DMSwarmSetCellDM(sw, dm, ierr);CHKERRA(ierr)
514cd81d59SMatthew G. Knepley      call DMSwarmRegisterPetscDatatypeField(sw, 'w_q', Nc, PETSC_SCALAR, ierr);CHKERRA(ierr)
524cd81d59SMatthew G. Knepley      call DMSwarmFinalizeFieldRegister(sw, ierr);CHKERRA(ierr)
534cd81d59SMatthew G. Knepley      call DMSwarmSetLocalSizes(sw, Np, zero, ierr);CHKERRA(ierr)
544cd81d59SMatthew G. Knepley      call DMSetFromOptions(sw, ierr);CHKERRA(ierr)
554cd81d59SMatthew G. Knepley      call DMSwarmGetField(sw, 'w_q', bs, dtype, wq, ierr);CHKERRA(ierr)
564cd81d59SMatthew G. Knepley      call DMSwarmGetField(sw, 'DMSwarmPIC_coor', bs, dtype, coords, ierr);CHKERRA(ierr)
574cd81d59SMatthew G. Knepley      do p = 1, Np
584cd81d59SMatthew G. Knepley        coords(p*2-1) = -cos(dble(p)/dble(Np+1) * PETSC_PI)
594cd81d59SMatthew G. Knepley        coords(p*2-0) =  sin(dble(p)/dble(Np+1) * PETSC_PI)
604cd81d59SMatthew G. Knepley        wq(p)         = 1.0
614cd81d59SMatthew G. Knepley      end do
624cd81d59SMatthew G. Knepley      call DMSwarmRestoreField(sw, 'DMSwarmPIC_coor', bs, dtype, coords, ierr);CHKERRA(ierr)
634cd81d59SMatthew G. Knepley      call DMSwarmRestoreField(sw, 'w_q', bs, dtype, wq, ierr);CHKERRA(ierr)
644cd81d59SMatthew G. Knepley      call DMSwarmMigrate(sw, removePoints, ierr);CHKERRA(ierr)
654cd81d59SMatthew G. Knepley      call DMViewFromOptions(sw, PETSC_NULL_VEC, '-swarm_view', ierr);CHKERRA(ierr)
664cd81d59SMatthew G. Knepley
674cd81d59SMatthew G. Knepley!     Project particles to field
684cd81d59SMatthew G. Knepley!       This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE
694cd81d59SMatthew G. Knepley      call DMCreateMassMatrix(sw, dm, M_p, ierr);CHKERRA(ierr)
704cd81d59SMatthew G. Knepley      call DMCreateGlobalVector(dm, rho, ierr);CHKERRA(ierr)
714cd81d59SMatthew G. Knepley      call DMSwarmCreateGlobalVectorFromField(sw, 'w_q', f, ierr);CHKERRA(ierr)
724cd81d59SMatthew G. Knepley      call MatMultTranspose(M_p, f, rho, ierr);CHKERRA(ierr)
734cd81d59SMatthew G. Knepley
744cd81d59SMatthew G. Knepley!     Visualize mesh field
754cd81d59SMatthew G. Knepley      call DMSetOutputSequenceNumber(dm, timestep, time, ierr);CHKERRA(ierr)
764cd81d59SMatthew G. Knepley      call PetscObjectViewFromOptions(rho, PETSC_NULL_VEC, '-rho_view', ierr);CHKERRA(ierr)
774cd81d59SMatthew G. Knepley
784cd81d59SMatthew G. Knepley!     Project field to particles
794cd81d59SMatthew G. Knepley!       This gives f_p = M_p^+ M f
804cd81d59SMatthew G. Knepley      call DMCreateMassMatrix(dm, dm, M, ierr);CHKERRA(ierr)
814cd81d59SMatthew G. Knepley      call DMCreateGlobalVector(dm, rhs, ierr);CHKERRA(ierr)
824cd81d59SMatthew G. Knepley      if (.false.) then
834cd81d59SMatthew G. Knepley         call MatMult(M, rho, rhs, ierr);CHKERRA(ierr) ! this is what you would do for and FE solve
844cd81d59SMatthew G. Knepley      else
854cd81d59SMatthew G. Knepley         call VecCopy(rho, rhs, ierr);CHKERRA(ierr) ! Indentity: M^1 M rho
864cd81d59SMatthew G. Knepley      end if
874cd81d59SMatthew G. Knepley      call KSPCreate(PETSC_COMM_WORLD, ksp, ierr);CHKERRA(ierr)
884cd81d59SMatthew G. Knepley      call KSPSetOptionsPrefix(ksp, 'ftop_', ierr);CHKERRA(ierr)
894cd81d59SMatthew G. Knepley      call KSPSetFromOptions(ksp, ierr);CHKERRA(ierr)
904cd81d59SMatthew G. Knepley      call KSPSetOperators(ksp, M_p, M_p, ierr);CHKERRA(ierr)
914cd81d59SMatthew G. Knepley      call KSPSolveTranspose(ksp, rhs, f, ierr);CHKERRA(ierr)
924cd81d59SMatthew G. Knepley      call KSPDestroy(ksp, ierr);CHKERRA(ierr)
934cd81d59SMatthew G. Knepley      call VecDestroy(rhs, ierr);CHKERRA(ierr)
944cd81d59SMatthew G. Knepley      call MatDestroy(M_p, ierr);CHKERRA(ierr)
954cd81d59SMatthew G. Knepley      call MatDestroy(M, ierr);CHKERRA(ierr)
964cd81d59SMatthew G. Knepley
974cd81d59SMatthew G. Knepley!     Visualize particle field
984cd81d59SMatthew G. Knepley      call DMSetOutputSequenceNumber(sw, timestep, time, ierr);CHKERRA(ierr)
994cd81d59SMatthew G. Knepley      call PetscObjectViewFromOptions(f, PETSC_NULL_VEC, '-weights_view', ierr);CHKERRA(ierr)
1004cd81d59SMatthew G. Knepley      call VecNorm(f,NORM_1,norm,ierr);CHKERRA(ierr)
1014cd81d59SMatthew G. Knepley      print *, 'Total number density = ', norm
1024cd81d59SMatthew G. Knepley!     Cleanup
1034cd81d59SMatthew G. Knepley      call DMSwarmDestroyGlobalVectorFromField(sw, 'w_q', f, ierr);CHKERRA(ierr)
1044cd81d59SMatthew G. Knepley      call MatDestroy(M_p, ierr);CHKERRA(ierr)
1054cd81d59SMatthew G. Knepley      call VecDestroy(rho, ierr);CHKERRA(ierr)
1064cd81d59SMatthew G. Knepley      call DMDestroy(sw, ierr);CHKERRA(ierr)
1074cd81d59SMatthew G. Knepley      call DMDestroy(dm, ierr);CHKERRA(ierr)
1084cd81d59SMatthew G. Knepley
1094cd81d59SMatthew G. Knepley      call PetscFinalize(ierr)
1104cd81d59SMatthew G. Knepley      end program DMSwarmTestProjection
1114cd81d59SMatthew G. Knepley
1124cd81d59SMatthew G. Knepley!/*TEST
1134cd81d59SMatthew G. Knepley!  build:
114*dfd57a17SPierre Jolivet!    requires: defined(PETSC_USING_F90FREEFORM) double !complex
1154cd81d59SMatthew G. Knepley!
1164cd81d59SMatthew G. Knepley!  test:
1174cd81d59SMatthew G. Knepley!    suffix: 0
1184cd81d59SMatthew G. Knepley!    requires: double
11930602db0SMatthew 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
1204cd81d59SMatthew G. Knepley!    filter: grep -v DM_ | grep -v atomic
1214cd81d59SMatthew G. Knepley!    filter_output: grep -v atomic
1224cd81d59SMatthew G. Knepley!
1234cd81d59SMatthew G. Knepley!TEST*/
124