1*4cd81d59SMatthew G. Knepley! Example program demonstrating projection between particle and finite element spaces 2*4cd81d59SMatthew G. Knepley program DMSwarmTestProjection 3*4cd81d59SMatthew G. Knepley#include "petsc/finclude/petscdmplex.h" 4*4cd81d59SMatthew G. Knepley#include "petsc/finclude/petscdmswarm.h" 5*4cd81d59SMatthew G. Knepley#include "petsc/finclude/petscksp.h" 6*4cd81d59SMatthew G. Knepley use petscdmplex 7*4cd81d59SMatthew G. Knepley use petscdmswarm 8*4cd81d59SMatthew G. Knepley use petscdt 9*4cd81d59SMatthew G. Knepley use petscksp 10*4cd81d59SMatthew G. Knepley use petscsys 11*4cd81d59SMatthew G. Knepley implicit none 12*4cd81d59SMatthew G. Knepley 13*4cd81d59SMatthew G. Knepley DM :: dm, sw 14*4cd81d59SMatthew G. Knepley PetscFE :: fe 15*4cd81d59SMatthew G. Knepley KSP :: ksp 16*4cd81d59SMatthew G. Knepley Mat :: M_p, M 17*4cd81d59SMatthew G. Knepley Vec :: f, rho, rhs 18*4cd81d59SMatthew G. Knepley PetscInt :: dim = 2, Nc = 1, degree = 1, timestep = 0 19*4cd81d59SMatthew G. Knepley PetscInt :: Np = 100, p, field = 0, zero = 0, bs 20*4cd81d59SMatthew G. Knepley PetscReal :: time = 0.0, norm 21*4cd81d59SMatthew G. Knepley PetscBool :: removePoints = PETSC_TRUE 22*4cd81d59SMatthew G. Knepley PetscDataType :: dtype 23*4cd81d59SMatthew G. Knepley PetscScalar, pointer :: coords(:) 24*4cd81d59SMatthew G. Knepley PetscScalar, pointer :: wq(:) 25*4cd81d59SMatthew G. Knepley PetscErrorCode :: ierr 26*4cd81d59SMatthew G. Knepley 27*4cd81d59SMatthew G. Knepley call PetscInitialize(PETSC_NULL_CHARACTER, ierr) 28*4cd81d59SMatthew G. Knepley if (ierr .ne. 0) then 29*4cd81d59SMatthew G. Knepley print*,'Unable to initialize PETSc' 30*4cd81d59SMatthew G. Knepley stop 31*4cd81d59SMatthew G. Knepley endif 32*4cd81d59SMatthew G. Knepley 33*4cd81d59SMatthew G. Knepley call DMPlexCreateBoxMesh(PETSC_COMM_WORLD, dim, PETSC_FALSE, PETSC_NULL_INTEGER, PETSC_NULL_REAL, PETSC_NULL_REAL, PETSC_NULL_INTEGER, PETSC_TRUE, dm, ierr);CHKERRA(ierr) 34*4cd81d59SMatthew G. Knepley call DMViewFromOptions(dm, PETSC_NULL_VEC, '-dm_view', ierr);CHKERRA(ierr) 35*4cd81d59SMatthew G. Knepley 36*4cd81d59SMatthew G. Knepley! Create finite element space 37*4cd81d59SMatthew G. Knepley call PetscFECreateLagrange(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, degree, PETSC_DETERMINE, fe, ierr);CHKERRA(ierr) 38*4cd81d59SMatthew G. Knepley call DMSetField(dm, field, PETSC_NULL_DMLABEL, fe, ierr);CHKERRA(ierr) 39*4cd81d59SMatthew G. Knepley call DMCreateDS(dm, ierr);CHKERRA(ierr) 40*4cd81d59SMatthew G. Knepley call PetscFEDestroy(fe, ierr);CHKERRA(ierr) 41*4cd81d59SMatthew G. Knepley 42*4cd81d59SMatthew G. Knepley! Create particle swarm 43*4cd81d59SMatthew G. Knepley call DMCreate(PETSC_COMM_WORLD, sw, ierr);CHKERRA(ierr) 44*4cd81d59SMatthew G. Knepley call DMSetType(sw, DMSWARM, ierr);CHKERRA(ierr) 45*4cd81d59SMatthew G. Knepley call DMSetDimension(sw, dim, ierr);CHKERRA(ierr) 46*4cd81d59SMatthew G. Knepley call DMSwarmSetType(sw, DMSWARM_PIC, ierr);CHKERRA(ierr) 47*4cd81d59SMatthew G. Knepley call DMSwarmSetCellDM(sw, dm, ierr);CHKERRA(ierr) 48*4cd81d59SMatthew G. Knepley call DMSwarmRegisterPetscDatatypeField(sw, 'w_q', Nc, PETSC_SCALAR, ierr);CHKERRA(ierr) 49*4cd81d59SMatthew G. Knepley call DMSwarmFinalizeFieldRegister(sw, ierr);CHKERRA(ierr) 50*4cd81d59SMatthew G. Knepley call DMSwarmSetLocalSizes(sw, Np, zero, ierr);CHKERRA(ierr) 51*4cd81d59SMatthew G. Knepley call DMSetFromOptions(sw, ierr);CHKERRA(ierr) 52*4cd81d59SMatthew G. Knepley call DMSwarmGetField(sw, 'w_q', bs, dtype, wq, ierr);CHKERRA(ierr) 53*4cd81d59SMatthew G. Knepley call DMSwarmGetField(sw, 'DMSwarmPIC_coor', bs, dtype, coords, ierr);CHKERRA(ierr) 54*4cd81d59SMatthew G. Knepley do p = 1, Np 55*4cd81d59SMatthew G. Knepley coords(p*2-1) = -cos(dble(p)/dble(Np+1) * PETSC_PI) 56*4cd81d59SMatthew G. Knepley coords(p*2-0) = sin(dble(p)/dble(Np+1) * PETSC_PI) 57*4cd81d59SMatthew G. Knepley wq(p) = 1.0 58*4cd81d59SMatthew G. Knepley end do 59*4cd81d59SMatthew G. Knepley call DMSwarmRestoreField(sw, 'DMSwarmPIC_coor', bs, dtype, coords, ierr);CHKERRA(ierr) 60*4cd81d59SMatthew G. Knepley call DMSwarmRestoreField(sw, 'w_q', bs, dtype, wq, ierr);CHKERRA(ierr) 61*4cd81d59SMatthew G. Knepley call DMSwarmMigrate(sw, removePoints, ierr);CHKERRA(ierr) 62*4cd81d59SMatthew G. Knepley call DMViewFromOptions(sw, PETSC_NULL_VEC, '-swarm_view', ierr);CHKERRA(ierr) 63*4cd81d59SMatthew G. Knepley 64*4cd81d59SMatthew G. Knepley! Project particles to field 65*4cd81d59SMatthew G. Knepley! This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE 66*4cd81d59SMatthew G. Knepley call DMCreateMassMatrix(sw, dm, M_p, ierr);CHKERRA(ierr) 67*4cd81d59SMatthew G. Knepley call DMCreateGlobalVector(dm, rho, ierr);CHKERRA(ierr) 68*4cd81d59SMatthew G. Knepley call DMSwarmCreateGlobalVectorFromField(sw, 'w_q', f, ierr);CHKERRA(ierr) 69*4cd81d59SMatthew G. Knepley call MatMultTranspose(M_p, f, rho, ierr);CHKERRA(ierr) 70*4cd81d59SMatthew G. Knepley 71*4cd81d59SMatthew G. Knepley! Visualize mesh field 72*4cd81d59SMatthew G. Knepley call DMSetOutputSequenceNumber(dm, timestep, time, ierr);CHKERRA(ierr) 73*4cd81d59SMatthew G. Knepley call PetscObjectViewFromOptions(rho, PETSC_NULL_VEC, '-rho_view', ierr);CHKERRA(ierr) 74*4cd81d59SMatthew G. Knepley 75*4cd81d59SMatthew G. Knepley! Project field to particles 76*4cd81d59SMatthew G. Knepley! This gives f_p = M_p^+ M f 77*4cd81d59SMatthew G. Knepley call DMCreateMassMatrix(dm, dm, M, ierr);CHKERRA(ierr) 78*4cd81d59SMatthew G. Knepley call DMCreateGlobalVector(dm, rhs, ierr);CHKERRA(ierr) 79*4cd81d59SMatthew G. Knepley if (.false.) then 80*4cd81d59SMatthew G. Knepley call MatMult(M, rho, rhs, ierr);CHKERRA(ierr) ! this is what you would do for and FE solve 81*4cd81d59SMatthew G. Knepley else 82*4cd81d59SMatthew G. Knepley call VecCopy(rho, rhs, ierr);CHKERRA(ierr) ! Indentity: M^1 M rho 83*4cd81d59SMatthew G. Knepley end if 84*4cd81d59SMatthew G. Knepley call KSPCreate(PETSC_COMM_WORLD, ksp, ierr);CHKERRA(ierr) 85*4cd81d59SMatthew G. Knepley call KSPSetOptionsPrefix(ksp, 'ftop_', ierr);CHKERRA(ierr) 86*4cd81d59SMatthew G. Knepley call KSPSetFromOptions(ksp, ierr);CHKERRA(ierr) 87*4cd81d59SMatthew G. Knepley call KSPSetOperators(ksp, M_p, M_p, ierr);CHKERRA(ierr) 88*4cd81d59SMatthew G. Knepley call KSPSolveTranspose(ksp, rhs, f, ierr);CHKERRA(ierr) 89*4cd81d59SMatthew G. Knepley call KSPDestroy(ksp, ierr);CHKERRA(ierr) 90*4cd81d59SMatthew G. Knepley call VecDestroy(rhs, ierr);CHKERRA(ierr) 91*4cd81d59SMatthew G. Knepley call MatDestroy(M_p, ierr);CHKERRA(ierr) 92*4cd81d59SMatthew G. Knepley call MatDestroy(M, ierr);CHKERRA(ierr) 93*4cd81d59SMatthew G. Knepley 94*4cd81d59SMatthew G. Knepley! Visualize particle field 95*4cd81d59SMatthew G. Knepley call DMSetOutputSequenceNumber(sw, timestep, time, ierr);CHKERRA(ierr) 96*4cd81d59SMatthew G. Knepley call PetscObjectViewFromOptions(f, PETSC_NULL_VEC, '-weights_view', ierr);CHKERRA(ierr) 97*4cd81d59SMatthew G. Knepley call VecNorm(f,NORM_1,norm,ierr);CHKERRA(ierr) 98*4cd81d59SMatthew G. Knepley print *, 'Total number density = ', norm 99*4cd81d59SMatthew G. Knepley! Cleanup 100*4cd81d59SMatthew G. Knepley call DMSwarmDestroyGlobalVectorFromField(sw, 'w_q', f, ierr);CHKERRA(ierr) 101*4cd81d59SMatthew G. Knepley call MatDestroy(M_p, ierr);CHKERRA(ierr) 102*4cd81d59SMatthew G. Knepley call VecDestroy(rho, ierr);CHKERRA(ierr) 103*4cd81d59SMatthew G. Knepley call DMDestroy(sw, ierr);CHKERRA(ierr) 104*4cd81d59SMatthew G. Knepley call DMDestroy(dm, ierr);CHKERRA(ierr) 105*4cd81d59SMatthew G. Knepley 106*4cd81d59SMatthew G. Knepley call PetscFinalize(ierr) 107*4cd81d59SMatthew G. Knepley end program DMSwarmTestProjection 108*4cd81d59SMatthew G. Knepley 109*4cd81d59SMatthew G. Knepley!/*TEST 110*4cd81d59SMatthew G. Knepley! build: 111*4cd81d59SMatthew G. Knepley! requires: define(PETSC_USING_F90FREEFORM) double !complex 112*4cd81d59SMatthew G. Knepley! 113*4cd81d59SMatthew G. Knepley! test: 114*4cd81d59SMatthew G. Knepley! suffix: 0 115*4cd81d59SMatthew G. Knepley! requires: double 116*4cd81d59SMatthew G. Knepley! args: -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 117*4cd81d59SMatthew G. Knepley! filter: grep -v DM_ | grep -v atomic 118*4cd81d59SMatthew G. Knepley! filter_output: grep -v atomic 119*4cd81d59SMatthew G. Knepley! 120*4cd81d59SMatthew G. Knepley!TEST*/ 121