13c611800SMark Adams static char help[] = "Example program demonstrating projection between particle and finite element spaces using OpenMP in 2D cylindrical coordinates\n"; 23c611800SMark Adams 33c611800SMark Adams #include "petscdmplex.h" 43c611800SMark Adams #include "petscds.h" 53c611800SMark Adams #include "petscdmswarm.h" 63c611800SMark Adams #include "petscksp.h" 73c611800SMark Adams #include <petsc/private/petscimpl.h> 83c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY) 93c611800SMark Adams #include <omp.h> 103c611800SMark Adams #endif 113c611800SMark Adams 123c611800SMark Adams typedef struct { 133c611800SMark Adams Mat MpTrans; 143c611800SMark Adams Mat Mp; 153c611800SMark Adams Vec ff; 163c611800SMark Adams Vec uu; 173c611800SMark Adams } MatShellCtx; 183c611800SMark Adams 193c611800SMark Adams PetscErrorCode MatMultMtM_SeqAIJ(Mat MtM,Vec xx,Vec yy) 203c611800SMark Adams { 213c611800SMark Adams MatShellCtx *matshellctx; 223c611800SMark Adams PetscErrorCode ierr; 233c611800SMark Adams 243c611800SMark Adams PetscFunctionBeginUser; 253c611800SMark Adams ierr = MatShellGetContext(MtM,&matshellctx);CHKERRQ(ierr); 26*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!matshellctx,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context"); 273c611800SMark Adams ierr = MatMult(matshellctx->Mp, xx, matshellctx->ff);CHKERRQ(ierr); 283c611800SMark Adams ierr = MatMult(matshellctx->MpTrans, matshellctx->ff, yy);CHKERRQ(ierr); 293c611800SMark Adams PetscFunctionReturn(0); 303c611800SMark Adams } 313c611800SMark Adams 323c611800SMark Adams PetscErrorCode MatMultAddMtM_SeqAIJ(Mat MtM,Vec xx, Vec yy, Vec zz) 333c611800SMark Adams { 343c611800SMark Adams MatShellCtx *matshellctx; 353c611800SMark Adams PetscErrorCode ierr; 363c611800SMark Adams 373c611800SMark Adams PetscFunctionBeginUser; 383c611800SMark Adams ierr = MatShellGetContext(MtM,&matshellctx);CHKERRQ(ierr); 39*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(!matshellctx,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No context"); 403c611800SMark Adams ierr = MatMult(matshellctx->Mp, xx, matshellctx->ff);CHKERRQ(ierr); 413c611800SMark Adams ierr = MatMultAdd(matshellctx->MpTrans, matshellctx->ff, yy, zz);CHKERRQ(ierr); 423c611800SMark Adams PetscFunctionReturn(0); 433c611800SMark Adams } 443c611800SMark Adams 453c611800SMark Adams PetscErrorCode createSwarm(const DM dm, DM *sw) 463c611800SMark Adams { 473c611800SMark Adams PetscErrorCode ierr; 483c611800SMark Adams PetscInt Nc = 1, dim = 2; 493c611800SMark Adams 503c611800SMark Adams PetscFunctionBeginUser; 513c611800SMark Adams ierr = DMCreate(PETSC_COMM_SELF, sw);CHKERRQ(ierr); 523c611800SMark Adams ierr = DMSetType(*sw, DMSWARM);CHKERRQ(ierr); 533c611800SMark Adams ierr = DMSetDimension(*sw, dim);CHKERRQ(ierr); 543c611800SMark Adams ierr = DMSwarmSetType(*sw, DMSWARM_PIC);CHKERRQ(ierr); 553c611800SMark Adams ierr = DMSwarmSetCellDM(*sw, dm);CHKERRQ(ierr); 563c611800SMark Adams ierr = DMSwarmRegisterPetscDatatypeField(*sw, "w_q", Nc, PETSC_SCALAR);CHKERRQ(ierr); 573c611800SMark Adams ierr = DMSwarmFinalizeFieldRegister(*sw);CHKERRQ(ierr); 583c611800SMark Adams ierr = DMSetFromOptions(*sw);CHKERRQ(ierr); 593c611800SMark Adams PetscFunctionReturn(0); 603c611800SMark Adams } 613c611800SMark Adams 623c611800SMark Adams PetscErrorCode gridToParticles(const DM dm, DM sw, PetscReal *moments, Vec rhs, Mat M_p) 633c611800SMark Adams { 643c611800SMark Adams PetscBool is_lsqr; 653c611800SMark Adams KSP ksp; 663c611800SMark Adams Mat PM_p=NULL,MtM,D; 673c611800SMark Adams Vec ff; 683c611800SMark Adams PetscErrorCode ierr; 693c611800SMark Adams PetscInt Np, timestep = 0, bs, N, M, nzl; 703c611800SMark Adams PetscReal time = 0.0; 713c611800SMark Adams PetscDataType dtype; 723c611800SMark Adams MatShellCtx *matshellctx; 733c611800SMark Adams 743c611800SMark Adams PetscFunctionBeginUser; 753c611800SMark Adams ierr = KSPCreate(PETSC_COMM_SELF, &ksp);CHKERRQ(ierr); 763c611800SMark Adams ierr = KSPSetOptionsPrefix(ksp, "ftop_");CHKERRQ(ierr); 773c611800SMark Adams ierr = KSPSetFromOptions(ksp);CHKERRQ(ierr); 783c611800SMark Adams ierr = PetscObjectTypeCompare((PetscObject)ksp,KSPLSQR,&is_lsqr); 793c611800SMark Adams if (!is_lsqr) { 803c611800SMark Adams ierr = MatGetLocalSize(M_p, &M, &N);CHKERRQ(ierr); 813c611800SMark Adams if (N>M) { 823c611800SMark Adams PC pc; 837d3de750SJacob Faibussowitsch ierr = PetscInfo(ksp, " M (%D) < M (%D) -- skip revert to lsqr\n",M,N);CHKERRQ(ierr); 843c611800SMark Adams is_lsqr = PETSC_TRUE; 853c611800SMark Adams ierr = KSPSetType(ksp,KSPLSQR);CHKERRQ(ierr); 863c611800SMark Adams ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 873c611800SMark Adams ierr = PCSetType(pc,PCNONE);CHKERRQ(ierr); // could put in better solver -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero 883c611800SMark Adams } else { 893c611800SMark Adams ierr = PetscNew(&matshellctx);CHKERRQ(ierr); 903c611800SMark Adams ierr = MatCreateShell(PetscObjectComm((PetscObject)dm),N,N,PETSC_DECIDE,PETSC_DECIDE,matshellctx,&MtM);CHKERRQ(ierr); 913c611800SMark Adams ierr = MatTranspose(M_p,MAT_INITIAL_MATRIX,&matshellctx->MpTrans);CHKERRQ(ierr); 923c611800SMark Adams matshellctx->Mp = M_p; 933c611800SMark Adams ierr = MatShellSetOperation(MtM, MATOP_MULT, (void (*)(void))MatMultMtM_SeqAIJ);CHKERRQ(ierr); 943c611800SMark Adams ierr = MatShellSetOperation(MtM, MATOP_MULT_ADD, (void (*)(void))MatMultAddMtM_SeqAIJ);CHKERRQ(ierr); 953c611800SMark Adams ierr = MatCreateVecs(M_p,&matshellctx->uu,&matshellctx->ff);CHKERRQ(ierr); 963c611800SMark Adams ierr = MatCreateSeqAIJ(PETSC_COMM_SELF,N,N,1,NULL,&D);CHKERRQ(ierr); 973c611800SMark Adams for (int i=0 ; i<N ; i++) { 983c611800SMark Adams const PetscScalar *vals; 993c611800SMark Adams const PetscInt *cols; 1003c611800SMark Adams PetscScalar dot = 0; 1013c611800SMark Adams ierr = MatGetRow(matshellctx->MpTrans,i,&nzl,&cols,&vals);CHKERRQ(ierr); 1023c611800SMark Adams for (int ii=0 ; ii<nzl ; ii++) dot += PetscSqr(vals[ii]); 103*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(dot==0.0,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Row %D is empty", i); 1043c611800SMark Adams ierr = MatSetValue(D,i,i,dot,INSERT_VALUES); 1053c611800SMark Adams } 1063c611800SMark Adams ierr = MatAssemblyBegin(D, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1073c611800SMark Adams ierr = MatAssemblyEnd(D, MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1087d3de750SJacob Faibussowitsch PetscInfo(M_p,"createMtMKSP Have %D eqs, nzl = %D\n",N,nzl); 1093c611800SMark Adams ierr = KSPSetOperators(ksp, MtM, D);CHKERRQ(ierr); 1103c611800SMark Adams ierr = MatViewFromOptions(D,NULL,"-ftop2_D_mat_view");CHKERRQ(ierr); 1113c611800SMark Adams ierr = MatViewFromOptions(M_p,NULL,"-ftop2_Mp_mat_view");CHKERRQ(ierr); 1123c611800SMark Adams ierr = MatViewFromOptions(matshellctx->MpTrans,NULL,"-ftop2_MpTranspose_mat_view");CHKERRQ(ierr); 1133c611800SMark Adams } 1143c611800SMark Adams } 1153c611800SMark Adams if (is_lsqr) { 1163c611800SMark Adams PC pc; 1173c611800SMark Adams PetscBool is_bjac; 1183c611800SMark Adams ierr = KSPGetPC(ksp,&pc);CHKERRQ(ierr); 1193c611800SMark Adams ierr = PetscObjectTypeCompare((PetscObject)pc,PCBJACOBI,&is_bjac); 1203c611800SMark Adams if (is_bjac) { 1213c611800SMark Adams ierr = DMSwarmCreateMassMatrixSquare(sw, dm, &PM_p);CHKERRQ(ierr); 1223c611800SMark Adams ierr = KSPSetOperators(ksp, M_p, PM_p);CHKERRQ(ierr); 1233c611800SMark Adams } else { 1243c611800SMark Adams ierr = KSPSetOperators(ksp, M_p, M_p);CHKERRQ(ierr); 1253c611800SMark Adams } 1263c611800SMark Adams } 1273c611800SMark Adams ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff);CHKERRQ(ierr); // this grabs access !!!!! 1283c611800SMark Adams if (!is_lsqr) { 1293c611800SMark Adams ierr = KSPSolve(ksp, rhs, matshellctx->uu);CHKERRQ(ierr); 1303c611800SMark Adams ierr = MatMult(M_p, matshellctx->uu, ff);CHKERRQ(ierr); 1313c611800SMark Adams ierr = MatDestroy(&matshellctx->MpTrans);CHKERRQ(ierr); 1323c611800SMark Adams ierr = VecDestroy(&matshellctx->ff);CHKERRQ(ierr); 1333c611800SMark Adams ierr = VecDestroy(&matshellctx->uu);CHKERRQ(ierr); 1343c611800SMark Adams ierr = MatDestroy(&D);CHKERRQ(ierr); 1353c611800SMark Adams ierr = MatDestroy(&MtM);CHKERRQ(ierr); 1363c611800SMark Adams ierr = PetscFree(matshellctx);CHKERRQ(ierr); 1373c611800SMark Adams } else { 1383c611800SMark Adams ierr = KSPSolveTranspose(ksp, rhs, ff);CHKERRQ(ierr); 1393c611800SMark Adams } 1403c611800SMark Adams ierr = KSPDestroy(&ksp);CHKERRQ(ierr); 1413c611800SMark Adams /* Visualize particle field */ 1423c611800SMark Adams ierr = DMSetOutputSequenceNumber(sw, timestep, time);CHKERRQ(ierr); 1433c611800SMark Adams ierr = VecViewFromOptions(ff, NULL, "-weights_view");CHKERRQ(ierr); 1443c611800SMark Adams ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff);CHKERRQ(ierr); 1453c611800SMark Adams 1463c611800SMark Adams /* compute energy */ 1473c611800SMark Adams if (moments) { 1483c611800SMark Adams PetscReal *wq, *coords; 1493c611800SMark Adams ierr = DMSwarmGetLocalSize(sw,&Np);CHKERRQ(ierr); 1503c611800SMark Adams ierr = DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr); 1513c611800SMark Adams ierr = DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr); 1523c611800SMark Adams moments[0] = moments[1] = moments[2] = 0; 1533c611800SMark Adams for (int p=0;p<Np;p++) { 1543c611800SMark Adams moments[0] += wq[p]; 1553c611800SMark Adams moments[1] += wq[p] * coords[p*2+0]; // x-momentum 1563c611800SMark Adams moments[2] += wq[p] * (PetscSqr(coords[p*2+0])+PetscSqr(coords[p*2+1])); 1573c611800SMark Adams } 1583c611800SMark Adams ierr = DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr); 1593c611800SMark Adams ierr = DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr); 1603c611800SMark Adams } 1613c611800SMark Adams ierr = MatDestroy(&PM_p); 1623c611800SMark Adams 1633c611800SMark Adams PetscFunctionReturn(0); 1643c611800SMark Adams } 1653c611800SMark Adams 1663c611800SMark Adams PetscErrorCode particlesToGrid(const DM dm, DM sw, const PetscInt Np, const PetscInt a_tid, const PetscInt dim, const PetscInt target, 1673c611800SMark Adams const PetscReal xx[], const PetscReal yy[], const PetscReal a_wp[], Vec rho, Mat *Mp_out) 1683c611800SMark Adams { 1693c611800SMark Adams 1703c611800SMark Adams PetscBool removePoints = PETSC_TRUE; 1713c611800SMark Adams PetscReal *wq, *coords; 1723c611800SMark Adams PetscDataType dtype; 1733c611800SMark Adams Mat M_p; 1743c611800SMark Adams Vec ff; 1753c611800SMark Adams PetscErrorCode ierr; 1763c611800SMark Adams PetscInt bs,p,zero=0; 1773c611800SMark Adams 1783c611800SMark Adams PetscFunctionBeginUser; 1793c611800SMark Adams ierr = DMSwarmSetLocalSizes(sw, Np, zero);CHKERRQ(ierr); 1803c611800SMark Adams ierr = DMSwarmGetField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr); 1813c611800SMark Adams ierr = DMSwarmGetField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr); 1823c611800SMark Adams for (p=0;p<Np;p++) { 1833c611800SMark Adams coords[p*2+0] = xx[p]; 1843c611800SMark Adams coords[p*2+1] = yy[p]; 1853c611800SMark Adams wq[p] = a_wp[p]; 1863c611800SMark Adams } 1873c611800SMark Adams ierr = DMSwarmRestoreField(sw, "DMSwarmPIC_coor", &bs, &dtype, (void**)&coords);CHKERRQ(ierr); 1883c611800SMark Adams ierr = DMSwarmRestoreField(sw, "w_q", &bs, &dtype, (void**)&wq);CHKERRQ(ierr); 1893c611800SMark Adams ierr = DMSwarmMigrate(sw, removePoints);CHKERRQ(ierr); 1903c611800SMark Adams ierr = PetscObjectSetName((PetscObject)sw, "Particle Grid");CHKERRQ(ierr); 1913c611800SMark Adams if (a_tid==target) {ierr = DMViewFromOptions(sw, NULL, "-swarm_view");CHKERRQ(ierr);} 1923c611800SMark Adams 1933c611800SMark Adams /* Project particles to field */ 1943c611800SMark Adams /* This gives M f = \int_\Omega \phi f, which looks like a rhs for a PDE */ 1953c611800SMark Adams ierr = DMCreateMassMatrix(sw, dm, &M_p);CHKERRQ(ierr); 1963c611800SMark Adams ierr = PetscObjectSetName((PetscObject)rho, "rho");CHKERRQ(ierr); 1973c611800SMark Adams 1983c611800SMark Adams ierr = DMSwarmCreateGlobalVectorFromField(sw, "w_q", &ff);CHKERRQ(ierr); // this grabs access !!!!! 1993c611800SMark Adams ierr = PetscObjectSetName((PetscObject)ff, "weights");CHKERRQ(ierr); 2003c611800SMark Adams ierr = MatMultTranspose(M_p, ff, rho);CHKERRQ(ierr); 2013c611800SMark Adams ierr = DMSwarmDestroyGlobalVectorFromField(sw, "w_q", &ff);CHKERRQ(ierr); 2023c611800SMark Adams 2033c611800SMark Adams /* Visualize mesh field */ 2043c611800SMark Adams if (a_tid==target) {ierr = VecViewFromOptions(rho, NULL, "-rho_view");CHKERRQ(ierr);} 2053c611800SMark Adams // output 2063c611800SMark Adams *Mp_out = M_p; 2073c611800SMark Adams 2083c611800SMark Adams PetscFunctionReturn(0); 2093c611800SMark Adams } 2103c611800SMark Adams static PetscErrorCode maxwellian(PetscInt dim, const PetscReal x[], PetscReal kt_m, PetscReal n, PetscScalar *u) 2113c611800SMark Adams { 2123c611800SMark Adams PetscInt i; 2133c611800SMark Adams PetscReal v2 = 0, theta = 2*kt_m; /* theta = 2kT/mc^2 */ 2143c611800SMark Adams 2153c611800SMark Adams PetscFunctionBegin; 2163c611800SMark Adams /* compute the exponents, v^2 */ 2173c611800SMark Adams for (i = 0; i < dim; ++i) v2 += x[i]*x[i]; 2183c611800SMark Adams /* evaluate the Maxwellian */ 2193c611800SMark Adams u[0] = n*PetscPowReal(PETSC_PI*theta,-1.5)*(PetscExpReal(-v2/theta)) * 2.*PETSC_PI*x[1]; // radial term for 2D axi-sym. 2203c611800SMark Adams PetscFunctionReturn(0); 2213c611800SMark Adams } 2223c611800SMark Adams #define NUM_SOLVE_LOOPS 100 2233c611800SMark Adams #define MAX_NUM_THRDS 12 2243c611800SMark Adams PetscErrorCode go() 2253c611800SMark Adams { 2263c611800SMark Adams DM dm_t[MAX_NUM_THRDS], sw_t[MAX_NUM_THRDS]; 2273c611800SMark Adams PetscFE fe; 2283c611800SMark Adams PetscInt dim = 2, Nc = 1, timestep = 0, i, faces[3]; 2293c611800SMark Adams PetscInt Np[2] = {10,10}, Np2[2], field = 0, target = 0, Np_t[MAX_NUM_THRDS]; 2303c611800SMark Adams PetscReal time = 0.0, moments_0[3], moments_1[3], vol; 2313c611800SMark Adams PetscReal lo[3] = {-5,0,-5}, hi[3] = {5,5,5}, h[3], hp[3], *xx_t[MAX_NUM_THRDS], *yy_t[MAX_NUM_THRDS], *wp_t[MAX_NUM_THRDS], solve_time = 0; 2323c611800SMark Adams Vec rho_t[MAX_NUM_THRDS], rhs_t[MAX_NUM_THRDS]; 2333c611800SMark Adams Mat M_p_t[MAX_NUM_THRDS]; 2343c611800SMark Adams PetscErrorCode ierr; 2353c611800SMark Adams #if defined PETSC_USE_LOG 2363c611800SMark Adams PetscLogStage stage; 2373c611800SMark Adams PetscLogEvent swarm_create_ev, solve_ev, solve_loop_ev; 2383c611800SMark Adams #endif 2393c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY) 2403c611800SMark Adams PetscInt numthreads = PetscNumOMPThreads; 2413c611800SMark Adams double starttime, endtime; 2423c611800SMark Adams #else 2433c611800SMark Adams PetscInt numthreads = 1; 2443c611800SMark Adams #endif 2453c611800SMark Adams 2463c611800SMark Adams PetscFunctionBeginUser; 2473c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY) 248*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(numthreads>MAX_NUM_THRDS,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "Too many threads %D > %D", numthreads, MAX_NUM_THRDS); 249*2c71b3e2SJacob Faibussowitsch PetscCheckFalse(numthreads<=0,PETSC_COMM_WORLD, PETSC_ERR_ARG_OUTOFRANGE, "No threads %D > %D ", numthreads, MAX_NUM_THRDS); 2503c611800SMark Adams #endif 2513c611800SMark Adams if (target >= numthreads) target = numthreads-1; 2523c611800SMark Adams ierr = PetscLogEventRegister("Create Swarm", DM_CLASSID, &swarm_create_ev);CHKERRQ(ierr); 2533c611800SMark Adams ierr = PetscLogEventRegister("Single solve", DM_CLASSID, &solve_ev);CHKERRQ(ierr); 2543c611800SMark Adams ierr = PetscLogEventRegister("Solve loop", DM_CLASSID, &solve_loop_ev);CHKERRQ(ierr); 2553c611800SMark Adams ierr = PetscLogStageRegister("Solve", &stage);CHKERRQ(ierr); 2563c611800SMark Adams i = dim; 2573c611800SMark Adams ierr = PetscOptionsGetIntArray(NULL, NULL, "-dm_plex_box_faces", faces, &i, NULL);CHKERRQ(ierr); 2583c611800SMark Adams i = dim; 2593c611800SMark Adams ierr = PetscOptionsGetIntArray(NULL, NULL, "-np", Np, &i, NULL);CHKERRQ(ierr); 2603c611800SMark Adams /* Create thread meshes */ 2613c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 2623c611800SMark Adams // setup mesh dm_t, could use PETSc's Landau create velocity space mesh here to get dm_t[tid] 2633c611800SMark Adams ierr = DMCreate(PETSC_COMM_SELF, &dm_t[tid]);CHKERRQ(ierr); 2643c611800SMark Adams ierr = DMSetType(dm_t[tid], DMPLEX);CHKERRQ(ierr); 2653c611800SMark Adams ierr = DMSetFromOptions(dm_t[tid]);CHKERRQ(ierr); 2663c611800SMark Adams ierr = PetscFECreateDefault(PETSC_COMM_SELF, dim, Nc, PETSC_FALSE, "", PETSC_DECIDE, &fe);CHKERRQ(ierr); 2673c611800SMark Adams ierr = PetscFESetFromOptions(fe);CHKERRQ(ierr); 2683c611800SMark Adams ierr = PetscObjectSetName((PetscObject)fe, "fe");CHKERRQ(ierr); 2693c611800SMark Adams ierr = DMSetField(dm_t[tid], field, NULL, (PetscObject)fe);CHKERRQ(ierr); 2703c611800SMark Adams ierr = DMCreateDS(dm_t[tid]);CHKERRQ(ierr); 2713c611800SMark Adams ierr = PetscFEDestroy(&fe);CHKERRQ(ierr); 2723c611800SMark Adams // helper vectors 2733c611800SMark Adams ierr = DMSetOutputSequenceNumber(dm_t[tid], timestep, time);CHKERRQ(ierr); // not used 2743c611800SMark Adams ierr = DMCreateGlobalVector(dm_t[tid], &rho_t[tid]);CHKERRQ(ierr); 2753c611800SMark Adams ierr = DMCreateGlobalVector(dm_t[tid], &rhs_t[tid]);CHKERRQ(ierr); 2763c611800SMark Adams // this mimics application code 2773c611800SMark Adams ierr = DMGetBoundingBox(dm_t[tid], lo, hi);CHKERRQ(ierr); 2783c611800SMark Adams if (tid==target) { 2793c611800SMark Adams ierr = DMViewFromOptions(dm_t[tid], NULL, "-dm_view");CHKERRQ(ierr); 2803c611800SMark Adams for (i=0,vol=1;i<dim;i++) { 2813c611800SMark Adams h[i] = (hi[i] - lo[i])/faces[i]; 2823c611800SMark Adams hp[i] = (hi[i] - lo[i])/Np[i]; 2833c611800SMark Adams vol *= (hi[i] - lo[i]); 2847d3de750SJacob Faibussowitsch ierr = PetscInfo(dm_t[tid]," lo = %g hi = %g n = %D h = %g hp = %g\n",lo[i],hi[i],faces[i],h[i],hp[i]);CHKERRQ(ierr); 2853c611800SMark Adams } 2863c611800SMark Adams } 2873c611800SMark Adams } 2883c611800SMark Adams // prepare particle data for problems. This mimics application code 2893c611800SMark Adams ierr = PetscLogEventBegin(swarm_create_ev,0,0,0,0);CHKERRQ(ierr); 2903c611800SMark Adams Np2[0] = Np[0]; Np2[1] = Np[1]; 2913c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { // change size of particle list a little 2923c611800SMark Adams Np_t[tid] = Np2[0]*Np2[1]; 2933c611800SMark Adams ierr = PetscMalloc3(Np_t[tid],&xx_t[tid],Np_t[tid],&yy_t[tid],Np_t[tid],&wp_t[tid]);CHKERRQ(ierr); 2943c611800SMark Adams if (tid==target) {moments_0[0] = moments_0[1] = moments_0[2] = 0;} 2953c611800SMark Adams for (int pi=0, pp=0;pi<Np2[0];pi++) { 2963c611800SMark Adams for (int pj=0;pj<Np2[1];pj++,pp++) { 2973c611800SMark Adams xx_t[tid][pp] = lo[0] + hp[0]/2. + pi*hp[0]; 2983c611800SMark Adams yy_t[tid][pp] = lo[1] + hp[1]/2. + pj*hp[1]; 2993c611800SMark Adams { 3003c611800SMark Adams PetscReal x[] = {xx_t[tid][pp],yy_t[tid][pp]}; 3013c611800SMark Adams ierr = maxwellian(2, x, 1.0, vol/(PetscReal)Np_t[tid], &wp_t[tid][pp]); 3023c611800SMark Adams } 3033c611800SMark Adams if (tid==target) { //energy_0 += wp_t[tid][pp]*(PetscSqr(xx_t[tid][pp])+PetscSqr(yy_t[tid][pp])); 3043c611800SMark Adams moments_0[0] += wp_t[tid][pp]; 3053c611800SMark Adams moments_0[1] += wp_t[tid][pp] * xx_t[tid][pp]; // x-momentum 3063c611800SMark Adams moments_0[2] += wp_t[tid][pp] * (PetscSqr(xx_t[tid][pp]) + PetscSqr(yy_t[tid][pp])); 3073c611800SMark Adams } 3083c611800SMark Adams } 3093c611800SMark Adams } 3103c611800SMark Adams Np2[0]++; Np2[1]++; 3113c611800SMark Adams } 3123c611800SMark Adams ierr = PetscLogEventEnd(swarm_create_ev,0,0,0,0);CHKERRQ(ierr); 3133c611800SMark Adams ierr = PetscLogEventBegin(solve_ev,0,0,0,0);CHKERRQ(ierr); 3143c611800SMark Adams /* Create particle swarm */ 3153c611800SMark Adams PetscPragmaOMP(parallel for) 3163c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3173c611800SMark Adams PetscErrorCode ierr_t; 3183c611800SMark Adams ierr_t = createSwarm(dm_t[tid], &sw_t[tid]); 3193c611800SMark Adams if (ierr_t) ierr = ierr_t; 3203c611800SMark Adams } 3213c611800SMark Adams CHKERRQ(ierr); 3223c611800SMark Adams PetscPragmaOMP(parallel for) 3233c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3243c611800SMark Adams PetscErrorCode ierr_t; 3253c611800SMark Adams ierr_t = particlesToGrid(dm_t[tid], sw_t[tid], Np_t[tid], tid, dim, target, xx_t[tid], yy_t[tid], wp_t[tid], rho_t[tid], &M_p_t[tid]); 3263c611800SMark Adams if (ierr_t) ierr = ierr_t; 3273c611800SMark Adams } 3283c611800SMark Adams CHKERRQ(ierr); 3293c611800SMark Adams /* Project field to particles */ 3303c611800SMark Adams /* This gives f_p = M_p^+ M f */ 3313c611800SMark Adams PetscPragmaOMP(parallel for) 3323c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3333c611800SMark Adams PetscErrorCode ierr_t; 3343c611800SMark Adams ierr_t = VecCopy(rho_t[tid], rhs_t[tid]); /* Identity: M^1 M rho */ 3353c611800SMark Adams if (ierr_t) ierr = ierr_t; 3363c611800SMark Adams } 3373c611800SMark Adams CHKERRQ(ierr); 3383c611800SMark Adams PetscPragmaOMP(parallel for) 3393c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3403c611800SMark Adams PetscErrorCode ierr_t; 3413c611800SMark Adams ierr_t = gridToParticles(dm_t[tid], sw_t[tid], (tid==target) ? moments_1 : NULL, rhs_t[tid], M_p_t[tid]); 3423c611800SMark Adams if (ierr_t) ierr = ierr_t; 3433c611800SMark Adams } 3443c611800SMark Adams CHKERRQ(ierr); 3453c611800SMark Adams /* Cleanup */ 3463c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3473c611800SMark Adams ierr = MatDestroy(&M_p_t[tid]);CHKERRQ(ierr); 3483c611800SMark Adams ierr = DMDestroy(&sw_t[tid]);CHKERRQ(ierr); 3493c611800SMark Adams } 3503c611800SMark Adams ierr = PetscLogEventEnd(solve_ev,0,0,0,0);CHKERRQ(ierr); 3513c611800SMark Adams /* for timing */ 3523c611800SMark Adams ierr = PetscOptionsClearValue(NULL,"-ftop_ksp_converged_reason");CHKERRQ(ierr); 3533c611800SMark Adams ierr = PetscOptionsClearValue(NULL,"-ftop_ksp_monitor");CHKERRQ(ierr); 3543c611800SMark Adams ierr = PetscOptionsClearValue(NULL,"-ftop_ksp_view");CHKERRQ(ierr); 3553c611800SMark Adams ierr = PetscOptionsClearValue(NULL,"-info");CHKERRQ(ierr); 3563c611800SMark Adams // repeat solve many times to get warmed up data 3573c611800SMark Adams ierr = PetscLogStagePush(stage);CHKERRQ(ierr); 3583c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY) 3593c611800SMark Adams starttime = MPI_Wtime(); 3603c611800SMark Adams #endif 3613c611800SMark Adams ierr = PetscLogEventBegin(solve_loop_ev,0,0,0,0);CHKERRQ(ierr); 3623c611800SMark Adams for (int d=0; d<NUM_SOLVE_LOOPS; d++) { 3633c611800SMark Adams /* Create particle swarm */ 3643c611800SMark Adams PetscPragmaOMP(parallel for) 3653c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3663c611800SMark Adams PetscErrorCode ierr_t; 3673c611800SMark Adams ierr_t = createSwarm(dm_t[tid], &sw_t[tid]); 3683c611800SMark Adams if (ierr_t) ierr = ierr_t; 3693c611800SMark Adams } 3703c611800SMark Adams CHKERRQ(ierr); 3713c611800SMark Adams PetscPragmaOMP(parallel for) 3723c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3733c611800SMark Adams PetscErrorCode ierr_t; 3743c611800SMark Adams ierr_t = particlesToGrid(dm_t[tid], sw_t[tid], Np_t[tid], tid, dim, target, xx_t[tid], yy_t[tid], wp_t[tid], rho_t[tid], &M_p_t[tid]); 3753c611800SMark Adams if (ierr_t) ierr = ierr_t; 3763c611800SMark Adams } 3773c611800SMark Adams CHKERRQ(ierr); 3783c611800SMark Adams PetscPragmaOMP(parallel for) 3793c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3803c611800SMark Adams PetscErrorCode ierr_t; 3813c611800SMark Adams ierr_t = VecCopy(rho_t[tid], rhs_t[tid]); /* Identity: M^1 M rho */ 3823c611800SMark Adams if (ierr_t) ierr = ierr_t; 3833c611800SMark Adams } 3843c611800SMark Adams CHKERRQ(ierr); 3853c611800SMark Adams PetscPragmaOMP(parallel for) 3863c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3873c611800SMark Adams PetscErrorCode ierr_t; 3883c611800SMark Adams ierr_t = gridToParticles(dm_t[tid], sw_t[tid], NULL, rhs_t[tid], M_p_t[tid]); 3893c611800SMark Adams if (ierr_t) ierr = ierr_t; 3903c611800SMark Adams } 3913c611800SMark Adams CHKERRQ(ierr); 3923c611800SMark Adams /* Cleanup */ 3933c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 3943c611800SMark Adams ierr = MatDestroy(&M_p_t[tid]);CHKERRQ(ierr); 3953c611800SMark Adams ierr = DMDestroy(&sw_t[tid]);CHKERRQ(ierr); 3963c611800SMark Adams } 3973c611800SMark Adams } 3983c611800SMark Adams ierr = PetscLogEventEnd(solve_loop_ev,0,0,0,0);CHKERRQ(ierr); 3993c611800SMark Adams #if defined(PETSC_HAVE_OPENMP) && defined(PETSC_HAVE_THREADSAFETY) 4003c611800SMark Adams endtime = MPI_Wtime(); 4013c611800SMark Adams solve_time += (endtime - starttime); 4023c611800SMark Adams #endif 4033c611800SMark Adams ierr = PetscLogStagePop();CHKERRQ(ierr); 4043c611800SMark Adams // 4053c611800SMark Adams ierr = PetscPrintf(PETSC_COMM_SELF,"Total number density: %20.12e (%20.12e); x-momentum = %g (%g); energy = %g error = %e, %D particles. Use %D threads, Solve time: %g\n", moments_1[0], moments_0[0], moments_1[1], moments_0[1], moments_1[2], (moments_1[2]-moments_0[2])/moments_0[2],Np[0]*Np[1],numthreads,solve_time);CHKERRQ(ierr); 4063c611800SMark Adams /* Cleanup */ 4073c611800SMark Adams for (int tid=0; tid<numthreads; tid++) { 4083c611800SMark Adams ierr = VecDestroy(&rho_t[tid]);CHKERRQ(ierr); 4093c611800SMark Adams ierr = VecDestroy(&rhs_t[tid]);CHKERRQ(ierr); 4103c611800SMark Adams ierr = DMDestroy(&dm_t[tid]);CHKERRQ(ierr); 4113c611800SMark Adams ierr = PetscFree3(xx_t[tid],yy_t[tid],wp_t[tid]);CHKERRQ(ierr); 4123c611800SMark Adams } 4133c611800SMark Adams 4143c611800SMark Adams PetscFunctionReturn(0); 4153c611800SMark Adams } 4163c611800SMark Adams 4173c611800SMark Adams int main(int argc, char **argv) 4183c611800SMark Adams { 4193c611800SMark Adams PetscErrorCode ierr; 4203c611800SMark Adams ierr = PetscInitialize(&argc, &argv, NULL,help);if (ierr) return ierr; 4213c611800SMark Adams ierr = go();CHKERRQ(ierr); 4223c611800SMark Adams ierr = PetscFinalize(); 4233c611800SMark Adams return ierr; 4243c611800SMark Adams } 4253c611800SMark Adams 4263c611800SMark Adams /*TEST 4273c611800SMark Adams 4283c611800SMark Adams build: 4293c611800SMark Adams requires: !complex 4303c611800SMark Adams 4313c611800SMark Adams test: 4323c611800SMark Adams suffix: 0 4333c611800SMark Adams requires: double triangle 4343c611800SMark Adams args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -petscspace_degree 2 -ftop_ksp_type lsqr -ftop_pc_type none -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14 4353c611800SMark Adams filter: grep -v DM_ | grep -v atomic 4363c611800SMark Adams 4373c611800SMark Adams test: 4383c611800SMark Adams suffix: 1 4393c611800SMark Adams requires: double triangle 4403c611800SMark Adams args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -petscspace_degree 2 -dm_plex_hash_location -ftop_ksp_type lsqr -ftop_pc_type bjacobi -ftop_sub_pc_type lu -ftop_sub_pc_factor_shift_type nonzero -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14 4413c611800SMark Adams filter: grep -v DM_ | grep -v atomic 4423c611800SMark Adams 4433c611800SMark Adams test: 4443c611800SMark Adams suffix: 2 4453c611800SMark Adams requires: double triangle 4463c611800SMark Adams args: -dm_plex_simplex 0 -dm_plex_box_faces 4,2 -dm_plex_box_lower -2.0,0.0 -dm_plex_box_upper 2.0,2.0 -petscspace_degree 2 -dm_plex_hash_location -ftop_ksp_type cg -ftop_pc_type jacobi -dm_view -ftop_ksp_converged_reason -ftop_ksp_rtol 1.e-14 4473c611800SMark Adams filter: grep -v DM_ | grep -v atomic 4483c611800SMark Adams 4493c611800SMark Adams TEST*/ 450