1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Tests DMSwarm\n\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown #include <petscdm.h> 5c4762a1bSJed Brown #include <petscdmda.h> 6c4762a1bSJed Brown #include <petscdmswarm.h> 7c4762a1bSJed Brown 8c4762a1bSJed Brown /* 9c4762a1bSJed Brown Checks for variable blocksize 10c4762a1bSJed Brown */ 11c4762a1bSJed Brown PetscErrorCode ex2_1(void) 12c4762a1bSJed Brown { 13c4762a1bSJed Brown DM dms; 14c4762a1bSJed Brown Vec x; 15c4762a1bSJed Brown PetscMPIInt rank; 16c4762a1bSJed Brown PetscInt p,bs,nlocal; 17c4762a1bSJed Brown 18362febeeSStefano Zampini PetscFunctionBegin; 19*9566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank)); 20c4762a1bSJed Brown 21*9566063dSJacob Faibussowitsch PetscCall(DMCreate(PETSC_COMM_WORLD,&dms)); 22*9566063dSJacob Faibussowitsch PetscCall(DMSetType(dms,DMSWARM)); 23*9566063dSJacob Faibussowitsch PetscCall(DMSwarmInitializeFieldRegister(dms)); 24*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"viscosity",1,PETSC_REAL)); 25*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRegisterPetscDatatypeField(dms,"strain",3,PETSC_REAL)); 26*9566063dSJacob Faibussowitsch PetscCall(DMSwarmFinalizeFieldRegister(dms)); 27*9566063dSJacob Faibussowitsch PetscCall(DMSwarmSetLocalSizes(dms,5+rank,4)); 28*9566063dSJacob Faibussowitsch PetscCall(DMView(dms,PETSC_VIEWER_STDOUT_WORLD)); 29*9566063dSJacob Faibussowitsch PetscCall(DMSwarmGetLocalSize(dms,&nlocal)); 30c4762a1bSJed Brown 31c4762a1bSJed Brown { 32c4762a1bSJed Brown PetscReal *array; 33*9566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"viscosity",&bs,NULL,(void**)&array)); 34c4762a1bSJed Brown for (p=0; p<nlocal; p++) { 35c4762a1bSJed Brown array[p] = 11.1 + p*0.1 + rank*100.0; 36c4762a1bSJed Brown } 37*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"viscosity",&bs,NULL,(void**)&array)); 38c4762a1bSJed Brown } 39c4762a1bSJed Brown 40c4762a1bSJed Brown { 41c4762a1bSJed Brown PetscReal *array; 42*9566063dSJacob Faibussowitsch PetscCall(DMSwarmGetField(dms,"strain",&bs,NULL,(void**)&array)); 43c4762a1bSJed Brown for (p=0; p<nlocal; p++) { 44c4762a1bSJed Brown array[bs*p+0] = 2.0e-2 + p*0.001 + rank*1.0; 45c4762a1bSJed Brown array[bs*p+1] = 2.0e-2 + p*0.002 + rank*1.0; 46c4762a1bSJed Brown array[bs*p+2] = 2.0e-2 + p*0.003 + rank*1.0; 47c4762a1bSJed Brown } 48*9566063dSJacob Faibussowitsch PetscCall(DMSwarmRestoreField(dms,"strain",&bs,NULL,(void**)&array)); 49c4762a1bSJed Brown } 50c4762a1bSJed Brown 51*9566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"viscosity",&x)); 52*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 53*9566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"viscosity",&x)); 54c4762a1bSJed Brown 55*9566063dSJacob Faibussowitsch PetscCall(DMSwarmCreateGlobalVectorFromField(dms,"strain",&x)); 56*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 57*9566063dSJacob Faibussowitsch PetscCall(DMSwarmDestroyGlobalVectorFromField(dms,"strain",&x)); 58c4762a1bSJed Brown 59*9566063dSJacob Faibussowitsch PetscCall(DMSwarmVectorDefineField(dms,"strain")); 60*9566063dSJacob Faibussowitsch PetscCall(DMCreateGlobalVector(dms,&x)); 61*9566063dSJacob Faibussowitsch PetscCall(VecView(x,PETSC_VIEWER_STDOUT_WORLD)); 62*9566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x)); 63*9566063dSJacob Faibussowitsch PetscCall(DMDestroy(&dms)); 64c4762a1bSJed Brown PetscFunctionReturn(0); 65c4762a1bSJed Brown } 66c4762a1bSJed Brown 67c4762a1bSJed Brown int main(int argc,char **argv) 68c4762a1bSJed Brown { 69c4762a1bSJed Brown 70*9566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc,&argv,(char*)0,help)); 71*9566063dSJacob Faibussowitsch PetscCall(ex2_1()); 72*9566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 73b122ec5aSJacob Faibussowitsch return 0; 74c4762a1bSJed Brown } 75c4762a1bSJed Brown 76c4762a1bSJed Brown /*TEST 77c4762a1bSJed Brown 78c4762a1bSJed Brown test: 79c4762a1bSJed Brown requires: !complex double 80c4762a1bSJed Brown nsize: 3 81c4762a1bSJed Brown filter: grep -v atomic 82c4762a1bSJed Brown filter_output: grep -v atomic 83c4762a1bSJed Brown 84c4762a1bSJed Brown TEST*/ 85