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 27*d8606c27SBarry Smith PetscCallA(PetscInitialize(PETSC_NULL_CHARACTER, ierr)) 28*d8606c27SBarry Smith PetscCallA(DMCreate(PETSC_COMM_WORLD, dm, ierr)) 29*d8606c27SBarry Smith PetscCallA(DMSetType(dm, DMPLEX, ierr)) 30*d8606c27SBarry Smith PetscCallA(DMSetFromOptions(dm, ierr)) 31*d8606c27SBarry Smith PetscCallA(DMGetDimension(dm, dim, ierr)) 32*d8606c27SBarry Smith PetscCallA(DMViewFromOptions(dm, PETSC_NULL_VEC, '-dm_view', ierr)) 334cd81d59SMatthew G. Knepley 344cd81d59SMatthew G. Knepley! Create finite element space 35*d8606c27SBarry Smith PetscCallA(PetscFECreateLagrange(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, degree, PETSC_DETERMINE, fe, ierr)) 36*d8606c27SBarry Smith PetscCallA(DMSetField(dm, field, PETSC_NULL_DMLABEL, fe, ierr)) 37*d8606c27SBarry Smith PetscCallA(DMCreateDS(dm, ierr)) 38*d8606c27SBarry Smith PetscCallA(PetscFEDestroy(fe, ierr)) 394cd81d59SMatthew G. Knepley 404cd81d59SMatthew G. Knepley! Create particle swarm 41*d8606c27SBarry Smith PetscCallA(DMCreate(PETSC_COMM_WORLD, sw, ierr)) 42*d8606c27SBarry Smith PetscCallA(DMSetType(sw, DMSWARM, ierr)) 43*d8606c27SBarry Smith PetscCallA(DMSetDimension(sw, dim, ierr)) 44*d8606c27SBarry Smith PetscCallA(DMSwarmSetType(sw, DMSWARM_PIC, ierr)) 45*d8606c27SBarry Smith PetscCallA(DMSwarmSetCellDM(sw, dm, ierr)) 46*d8606c27SBarry Smith PetscCallA(DMSwarmRegisterPetscDatatypeField(sw, 'w_q', Nc, PETSC_SCALAR, ierr)) 47*d8606c27SBarry Smith PetscCallA(DMSwarmFinalizeFieldRegister(sw, ierr)) 48*d8606c27SBarry Smith PetscCallA(DMSwarmSetLocalSizes(sw, Np, zero, ierr)) 49*d8606c27SBarry Smith PetscCallA(DMSetFromOptions(sw, ierr)) 50*d8606c27SBarry Smith PetscCallA(DMSwarmGetField(sw, 'w_q', bs, dtype, wq, ierr)) 51*d8606c27SBarry Smith PetscCallA(DMSwarmGetField(sw, 'DMSwarmPIC_coor', bs, dtype, coords, ierr)) 524cd81d59SMatthew G. Knepley do p = 1, Np 534cd81d59SMatthew G. Knepley coords(p*2-1) = -cos(dble(p)/dble(Np+1) * PETSC_PI) 544cd81d59SMatthew G. Knepley coords(p*2-0) = sin(dble(p)/dble(Np+1) * PETSC_PI) 554cd81d59SMatthew G. Knepley wq(p) = 1.0 564cd81d59SMatthew G. Knepley end do 57*d8606c27SBarry Smith PetscCallA(DMSwarmRestoreField(sw, 'DMSwarmPIC_coor', bs, dtype, coords, ierr)) 58*d8606c27SBarry Smith PetscCallA(DMSwarmRestoreField(sw, 'w_q', bs, dtype, wq, ierr)) 59*d8606c27SBarry Smith PetscCallA(DMSwarmMigrate(sw, removePoints, ierr)) 60*d8606c27SBarry Smith PetscCallA(DMViewFromOptions(sw, PETSC_NULL_VEC, '-swarm_view', ierr)) 614cd81d59SMatthew G. Knepley 624cd81d59SMatthew G. Knepley! Project particles to field 634cd81d59SMatthew G. Knepley! This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE 64*d8606c27SBarry Smith PetscCallA(DMCreateMassMatrix(sw, dm, M_p, ierr)) 65*d8606c27SBarry Smith PetscCallA(DMCreateGlobalVector(dm, rho, ierr)) 66*d8606c27SBarry Smith PetscCallA(DMSwarmCreateGlobalVectorFromField(sw, 'w_q', f, ierr)) 67*d8606c27SBarry Smith PetscCallA(MatMultTranspose(M_p, f, rho, ierr)) 684cd81d59SMatthew G. Knepley 694cd81d59SMatthew G. Knepley! Visualize mesh field 70*d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(dm, timestep, time, ierr)) 71*d8606c27SBarry Smith PetscCallA(PetscObjectViewFromOptions(rho, PETSC_NULL_VEC, '-rho_view', ierr)) 724cd81d59SMatthew G. Knepley 734cd81d59SMatthew G. Knepley! Project field to particles 744cd81d59SMatthew G. Knepley! This gives f_p = M_p^+ M f 75*d8606c27SBarry Smith PetscCallA(DMCreateMassMatrix(dm, dm, M, ierr)) 76*d8606c27SBarry Smith PetscCallA(DMCreateGlobalVector(dm, rhs, ierr)) 774cd81d59SMatthew G. Knepley if (.false.) then 78*d8606c27SBarry Smith PetscCallA(MatMult(M, rho, rhs, ierr)) ! this is what you would do for and FE solve 794cd81d59SMatthew G. Knepley else 80*d8606c27SBarry Smith PetscCallA(VecCopy(rho, rhs, ierr)) ! Identity: M^1 M rho 814cd81d59SMatthew G. Knepley end if 82*d8606c27SBarry Smith PetscCallA(KSPCreate(PETSC_COMM_WORLD, ksp, ierr)) 83*d8606c27SBarry Smith PetscCallA(KSPSetOptionsPrefix(ksp, 'ftop_', ierr)) 84*d8606c27SBarry Smith PetscCallA(KSPSetFromOptions(ksp, ierr)) 85*d8606c27SBarry Smith PetscCallA(KSPSetOperators(ksp, M_p, M_p, ierr)) 86*d8606c27SBarry Smith PetscCallA(KSPSolveTranspose(ksp, rhs, f, ierr)) 87*d8606c27SBarry Smith PetscCallA(KSPDestroy(ksp, ierr)) 88*d8606c27SBarry Smith PetscCallA(VecDestroy(rhs, ierr)) 89*d8606c27SBarry Smith PetscCallA(MatDestroy(M_p, ierr)) 90*d8606c27SBarry Smith PetscCallA(MatDestroy(M, ierr)) 914cd81d59SMatthew G. Knepley 924cd81d59SMatthew G. Knepley! Visualize particle field 93*d8606c27SBarry Smith PetscCallA(DMSetOutputSequenceNumber(sw, timestep, time, ierr)) 94*d8606c27SBarry Smith PetscCallA(PetscObjectViewFromOptions(f, PETSC_NULL_VEC, '-weights_view', ierr)) 95*d8606c27SBarry Smith PetscCallA(VecNorm(f,NORM_1,norm,ierr)) 964cd81d59SMatthew G. Knepley print *, 'Total number density = ', norm 974cd81d59SMatthew G. Knepley! Cleanup 98*d8606c27SBarry Smith PetscCallA(DMSwarmDestroyGlobalVectorFromField(sw, 'w_q', f, ierr)) 99*d8606c27SBarry Smith PetscCallA(MatDestroy(M_p, ierr)) 100*d8606c27SBarry Smith PetscCallA(VecDestroy(rho, ierr)) 101*d8606c27SBarry Smith PetscCallA(DMDestroy(sw, ierr)) 102*d8606c27SBarry Smith PetscCallA(DMDestroy(dm, ierr)) 1034cd81d59SMatthew G. Knepley 104*d8606c27SBarry Smith PetscCallA(PetscFinalize(ierr)) 1054cd81d59SMatthew G. Knepley end program DMSwarmTestProjection 1064cd81d59SMatthew G. Knepley 1074cd81d59SMatthew G. Knepley!/*TEST 1084cd81d59SMatthew G. Knepley! build: 109dfd57a17SPierre Jolivet! requires: defined(PETSC_USING_F90FREEFORM) double !complex 1104cd81d59SMatthew G. Knepley! 1114cd81d59SMatthew G. Knepley! test: 1124cd81d59SMatthew G. Knepley! suffix: 0 1134cd81d59SMatthew G. Knepley! requires: double 11430602db0SMatthew 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 1154cd81d59SMatthew G. Knepley! filter: grep -v DM_ | grep -v atomic 1164cd81d59SMatthew G. Knepley! filter_output: grep -v atomic 1174cd81d59SMatthew G. Knepley! 1184cd81d59SMatthew G. Knepley!TEST*/ 119