movedata.c

Go to the documentation of this file.
00001 /*
00002     This file is part of magpar.
00003 
00004     Copyright (C) 2002-2010 Werner Scholz
00005 
00006     www:   http://www.magpar.net/
00007     email: magpar(at)magpar.net
00008 
00009     magpar is free software; you can redistribute it and/or modify
00010     it under the terms of the GNU General Public License as published by
00011     the Free Software Foundation; either version 2 of the License, or
00012     (at your option) any later version.
00013 
00014     magpar is distributed in the hope that it will be useful,
00015     but WITHOUT ANY WARRANTY; without even the implied warranty of
00016     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
00017     GNU General Public License for more details.
00018 
00019     You should have received a copy of the GNU General Public License
00020     along with magpar; if not, write to the Free Software
00021     Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA  02111-1307  USA
00022 */
00023 
00024 static char const Id[] = "$Id: movedata.c 3054 2010-04-15 17:03:46Z scholz $\n\n";
00025 static char const Td[] = "$Today: " __FILE__ "  " __DATE__ "  "  __TIME__  " $\n\n";
00026 
00027 #include "init.h"
00028 #include "util/util.h"
00029 
00030 int DataMoveData(GridData *gdata)
00031 {
00032   MagparFunctionLogBegin;
00033 
00034   int rank,size;
00035   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
00036   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
00037 
00038   if (size==1) {
00039 #ifdef HAVE_ELEVERTALL
00040     gdata->elevertall=gdata->elevert;
00041 #endif
00042 
00043     int *ia,*ja;
00044 
00045     ierr = Mesh2Nodal(gdata->n_ele,gdata->n_vert,gdata->elevert,&ia,&ja);CHKERRQ(ierr);
00046     gdata->mesh2nodal_ia=ia;
00047     gdata->mesh2nodal_ja=ja;
00048 
00049     ierr = Mesh2Dual(gdata->n_ele,gdata->n_vert,gdata->elevert,&ia,&ja);CHKERRQ(ierr);
00050     gdata->mesh2dual_ia=ia;
00051     gdata->mesh2dual_ja=ja;
00052 
00053     MagparFunctionLogReturn(0);
00054   }
00055 
00056   /* broadcast vertxyz to all processors */
00057   ierr = MPI_Bcast(&gdata->n_vert,    1,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00058   if (rank) {
00059     ierr = PetscMalloc(ND*gdata->n_vert*sizeof(PetscReal),&gdata->vertxyz);CHKERRQ(ierr);
00060   }
00061   ierr = MPI_Bcast(gdata->vertxyz,gdata->n_vert*ND,MPIU_SCALAR,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00062 
00063   ierr = MPI_Bcast(&gdata->n_prop,1,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00064   if (rank) {
00065     ierr = PetscMalloc(gdata->n_prop*NP*sizeof(PetscReal),&gdata->propdat);CHKERRQ(ierr);
00066   }
00067   ierr = MPI_Bcast(gdata->propdat,gdata->n_prop*NP,MPIU_SCALAR,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00068 
00069   ierr = MPI_Bcast(&gdata->n_ele,1,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00070 
00071 #ifdef HAVE_ELEVERTALL
00072   ierr = PetscMalloc(NV*gdata->n_ele*sizeof(int),&gdata->elevertall);CHKERRQ(ierr);
00073   if (!rank) {
00074     ierr = PetscMemcpy(gdata->elevertall,gdata->elevert,NV*gdata->n_ele*sizeof(int));CHKERRQ(ierr)
00075   }
00076   ierr = MPI_Bcast(gdata->elevertall,NV*gdata->n_ele,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00077 #ifdef HAVE_M2IJ
00078   int *ia,*ja;
00079   if (!rank) {
00080     ierr = Mesh2Nodal(gdata->n_ele,gdata->n_vert,gdata->elevert,&ia,&ja);CHKERRQ(ierr);
00081     gdata->mesh2nodal_ia=ia;
00082     gdata->mesh2nodal_ja=ja;
00083 
00084     ierr = Mesh2Dual(gdata->n_ele,gdata->n_vert,gdata->elevert,&ia,&ja);CHKERRQ(ierr);
00085     gdata->mesh2dual_ia=ia;
00086     gdata->mesh2dual_ja=ja;
00087   }
00088   else {
00089     ierr = PetscMalloc((gdata->n_vert+1)*sizeof(int),&gdata->mesh2nodal_ia);CHKERRQ(ierr);
00090     ierr = PetscMalloc((gdata->n_ele+1)*sizeof(int),&gdata->mesh2dual_ia);CHKERRQ(ierr);
00091   }
00092   ierr = MPI_Bcast(gdata->mesh2nodal_ia,gdata->n_vert+1,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00093   ierr = MPI_Bcast(gdata->mesh2dual_ia,gdata->n_ele+1,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00094 
00095   if (rank) {
00096     ierr = PetscMalloc(gdata->mesh2nodal_ia[gdata->n_vert]*sizeof(int),&gdata->mesh2nodal_ja);CHKERRQ(ierr);
00097     ierr = PetscMalloc(gdata->mesh2dual_ia[gdata->n_ele]*sizeof(int),&gdata->mesh2dual_ja);CHKERRQ(ierr);
00098   }
00099   ierr = MPI_Bcast(gdata->mesh2nodal_ja,gdata->mesh2nodal_ia[gdata->n_vert],MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00100   ierr = MPI_Bcast(gdata->mesh2dual_ja,gdata->mesh2dual_ia[gdata->n_ele],MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00101 #endif
00102 #endif
00103 
00104   int t_gotele;
00105   t_gotele=ascat(&gdata->elevert,gdata->ln_ele,NV,gdata->elenewproc,size,1);
00106   int t_gotn;
00107   t_gotn=  ascat(&gdata->eleprop,gdata->ln_ele,1, gdata->elenewproc,size,1);assert(t_gotele==t_gotn);
00108 
00109   /* distribute boundary faces evenly over processors */
00110   ierr = MPI_Bcast(&gdata->n_vert_bnd,1,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00111   ierr = MPI_Bcast(&gdata->n_bnd_fac, 1,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00112   int *bndfacnewproc;
00113   ierr = PetscMalloc(gdata->n_bnd_fac*sizeof(int),&bndfacnewproc);CHKERRQ(ierr);
00114   if (!rank) {
00115     for (int i=0;i<gdata->n_bnd_fac;i++) {
00116       bndfacnewproc[i]=distint(gdata->n_bnd_fac,size,i);
00117     }
00118   }
00119   /* cannot move to facnb.c: partsurfser.c needs bndfacvert on proc 0 */
00120   int t_gotbndfac;
00121   t_gotbndfac=ascat(&gdata->bndfacvert,gdata->ln_bnd_fac,NN,bndfacnewproc,size,1);
00122   gdata->ln_bnd_fac=t_gotbndfac;
00123   ierr = PetscFree(bndfacnewproc);CHKERRQ(ierr);
00124   PetscSynchronizedPrintf(PETSC_COMM_WORLD,
00125     "<%i>Received local number of boundary faces: %i\n",
00126     rank,gdata->ln_bnd_fac
00127   );
00128   PetscSynchronizedFlush(PETSC_COMM_WORLD);
00129   int tsum;
00130   ierr=MPI_Allreduce((int*)&gdata->ln_bnd_fac,(int*)&tsum,1,MPI_INT,MPI_SUM,PETSC_COMM_WORLD);CHKERRQ(ierr);
00131   if (!rank) assert(gdata->n_bnd_fac==tsum);
00132   PetscPrintf(PETSC_COMM_WORLD,
00133     "Total number of boundary faces: %i\n",
00134     gdata->n_bnd_fac
00135   );
00136 
00137   /* Create a global reordering of the element numbers */
00138   int t_gotvert;
00139   t_gotvert=ascat(&gdata->vertl2g,gdata->ln_vert,1,gdata->vertnewproc,size,1);
00140   t_gotn=ascat(&gdata->elel2g,gdata->ln_ele,1,gdata->elenewproc,size,1);assert(t_gotele==t_gotn);
00141 
00142   /* broadcast all boundary vertex indicators to all processors
00143      (required for bmatrix)
00144   */
00145   if (rank) {
00146     ierr = PetscMalloc(gdata->n_vert*sizeof(int),&gdata->vertbndg2bnd);CHKERRQ(ierr);
00147   }
00148   ierr = MPI_Bcast(gdata->vertbndg2bnd,gdata->n_vert,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00149 
00150   t_gotn=ascat(&gdata->vertprop,gdata->ln_vert,1,gdata->vertnewproc,size,1);assert(t_gotvert==t_gotn);
00151 
00152   /* scatter magnetization from the first processor
00153      into the distributed vector
00154   */
00155 
00156   Vec tv_M;
00157   ierr = VecCreate(PETSC_COMM_WORLD,&tv_M);CHKERRQ(ierr);
00158   ierr = VecSetSizes(tv_M,ND*t_gotvert,ND*gdata->n_vert);CHKERRQ(ierr);
00159   ierr = VecSetFromOptions(tv_M);CHKERRQ(ierr);
00160   ierr = VecSetBlockSize(tv_M,ND);CHKERRQ(ierr);
00161 
00162   IS tis_from,tis_to;
00163   if (!rank) {
00164     /* preserve the ordering since the data have already been permuted above */
00165     ierr = ISCreateStride(PETSC_COMM_SELF,ND*gdata->n_vert,0,1,&tis_from);CHKERRQ(ierr);
00166     ierr = ISCreateStride(PETSC_COMM_SELF,ND*gdata->n_vert,0,1,&tis_to);CHKERRQ(ierr);
00167 
00168   }
00169   else {
00170     PetscInt t_idx_from[1],t_idx_to[1];
00171     tis_from=NULL;
00172     tis_to=NULL;
00173     ISCreateGeneral(PETSC_COMM_SELF,0,t_idx_from,&tis_from);
00174     ISCreateGeneral(PETSC_COMM_SELF,0,t_idx_to,  &tis_to);
00175   }
00176 
00177   VecScatter t_scatter;
00178   VecScatterCreate(gdata->M,tis_from,tv_M,tis_to,&t_scatter);
00179   ISDestroy(tis_from);
00180   ISDestroy(tis_to);
00181 
00182   VecScatterBegin(t_scatter,gdata->M,tv_M,INSERT_VALUES,SCATTER_FORWARD);
00183   VecScatterEnd  (t_scatter,gdata->M,tv_M,INSERT_VALUES,SCATTER_FORWARD);
00184   VecDestroy(gdata->M);
00185   gdata->M=tv_M;
00186 
00187 #ifdef ADDONS
00188   /* broadcast vertex coupling information to all processors */
00189   int j;
00190   if (!rank) {
00191     if (gdata->vertcoupl!=PETSC_NULL) {
00192       j=1;
00193     }
00194     else {
00195       j=0;
00196     }
00197   }
00198   ierr = MPI_Bcast(&j,1,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00199 
00200   if (j) {
00201     if (rank) {
00202       ierr = PetscMalloc(gdata->n_vert*sizeof(int),&gdata->vertcoupl);CHKERRQ(ierr);
00203     }
00204 
00205     ierr = MPI_Bcast(gdata->vertcoupl,gdata->n_vert,MPI_INT,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
00206   }
00207   else {
00208     gdata->vertcoupl=PETSC_NULL;
00209   }
00210 
00211   if (gdata->VH1!=PETSC_NULL) {
00212     Vec tv_H1;
00213     ierr = VecDuplicate(tv_M,&tv_H1);CHKERRQ(ierr);
00214     VecScatterBegin(t_scatter,gdata->VH1,tv_H1,INSERT_VALUES,SCATTER_FORWARD);
00215     VecScatterEnd  (t_scatter,gdata->VH1,tv_H1,INSERT_VALUES,SCATTER_FORWARD);
00216     VecDestroy(gdata->VH1);
00217     gdata->VH1=tv_H1;
00218   }
00219   if (gdata->VH2!=PETSC_NULL) {
00220     Vec tv_H2;
00221     ierr = VecDuplicate(tv_M,&tv_H2);CHKERRQ(ierr);
00222     VecScatterBegin(t_scatter,gdata->VH2,tv_H2,INSERT_VALUES,SCATTER_FORWARD);
00223     VecScatterEnd  (t_scatter,gdata->VH2,tv_H2,INSERT_VALUES,SCATTER_FORWARD);
00224     VecDestroy(gdata->VH2);
00225     gdata->VH2=tv_H2;
00226   }
00227 #endif
00228   VecScatterDestroy(t_scatter);
00229 
00230   t_gotn=ascat(&gdata->elenewproc,gdata->ln_ele,1,gdata->elenewproc,size,1);assert(t_gotele==t_gotn);
00231   t_gotn=ascat(&gdata->vertnewproc,gdata->ln_vert,1,gdata->vertnewproc,size,1);assert(t_gotvert==t_gotn);
00232 
00233   gdata->ln_ele=t_gotele;
00234   gdata->ln_vert=t_gotvert;
00235 
00236   PetscSynchronizedPrintf(PETSC_COMM_WORLD,"<%i>Received local number of elements: %i\n",rank,gdata->ln_ele);
00237   PetscSynchronizedPrintf(PETSC_COMM_WORLD,"<%i>Received local number of vertices: %i\n",rank,gdata->ln_vert);
00238   PetscSynchronizedFlush(PETSC_COMM_WORLD);
00239 
00240   /* check correct data transfer */
00241   for (int i=0; i<gdata->ln_ele; i++)
00242   {
00243     assert(gdata->elenewproc[i]==rank);
00244   }
00245   for (int i=0; i<gdata->ln_vert; i++)
00246   {
00247     assert(gdata->vertnewproc[i]==rank);
00248   }
00249   PetscSynchronizedFlush(PETSC_COMM_WORLD);
00250 
00251   int vidlow, vidhigh, vlocsize;
00252   ierr = VecGetOwnershipRange(gdata->M,(PetscInt*)&vidlow,(PetscInt*)&vidhigh);CHKERRQ(ierr);
00253   ierr = VecGetLocalSize(gdata->M,(PetscInt*)&vlocsize);CHKERRQ(ierr);
00254   vidlow   /= ND;
00255   vidhigh  /= ND;
00256   vlocsize /= ND;
00257   assert(vidhigh-vidlow==vlocsize);
00258   assert(vlocsize==gdata->ln_vert);
00259   PetscSynchronizedPrintf(PETSC_COMM_WORLD,"<%i>Ownership(M): from vertex id %i to vertex id %i\n",rank,vidlow,vidhigh-1);
00260   PetscSynchronizedFlush(PETSC_COMM_WORLD);
00261 
00262   ierr = PetscFree(gdata->vertnewproc);CHKERRQ(ierr);
00263   gdata->vertnewproc=PETSC_NULL;
00264   ierr = PetscFree(gdata->elenewproc);CHKERRQ(ierr);
00265   gdata->elenewproc=PETSC_NULL;
00266 
00267   MagparFunctionLogReturn(0);
00268 }
00269 

magpar - Parallel Finite Element Micromagnetics Package
Copyright (C) 2002-2010 Werner Scholz