filterelements.c

Go to the documentation of this file.
00001 /*
00002     This file is part of magpar.
00003 
00004     Copyright (C) 2002-2009 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: filterelements.c 2681 2009-07-31 04:30:53Z scholz $\n\n";
00025 static char const Td[] = "$Today: " __FILE__ "  " __DATE__ "  "  __TIME__  " $\n\n";
00026 
00027 #include "init.h"
00028 
00029 #define C_DEL -5   
00031 int FilterElements(int *pn_ele,int **pelevert,int **peleprop,PetscReal *propdat)
00032 {
00033   MagparFunctionLogBegin;
00034 
00035   int rank,size;
00036   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
00037   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
00038 
00039   if (rank) {
00040     MagparFunctionLogReturn(0);
00041   }
00042 
00043   int n_ele;    n_ele=*pn_ele;
00044   int *elevert; elevert=*pelevert;
00045   int *eleprop; eleprop=*peleprop;
00046 
00047   int n_del=0;
00048 
00049   /* get number of elements to be deleted, delete element if Js<0 */
00050   for (int i=0;i<n_ele;i++) {
00051     if (propdat[NP*eleprop[i]+4]<0.0) {
00052       eleprop[i]=-1;
00053       n_del++;
00054     }
00055   }
00056 
00057   if (n_del>0) {
00058     PetscPrintf(PETSC_COMM_WORLD,
00059     "Removing elements: %i out of %i\n",n_del,n_ele);
00060   }
00061   else {
00062     PetscPrintf(PETSC_COMM_WORLD,"No elements to be removed\n");
00063     MagparFunctionLogReturn(0);
00064   }
00065 
00066 
00067   PetscPrintf(PETSC_COMM_WORLD,"Allocating new element arrays\n");
00068 
00069   /* allocate new array for elements */
00070   int *televert;
00071   ierr = PetscMalloc(NV*(n_ele-n_del)*sizeof(int),&televert);CHKERRQ(ierr);
00072   int *teleprop;
00073   ierr = PetscMalloc((n_ele-n_del)*sizeof(int),&teleprop);CHKERRQ(ierr);
00074 
00075   PetscPrintf(PETSC_COMM_WORLD,"Copying element data\n");
00076 
00077   /* copy elements and properties */
00078   int k=0;
00079   for (int i=0;i<n_ele;i++) {
00080     if (eleprop[i]>=0) {
00081       for (int j=0;j<NV;j++) {
00082         televert[NV*k+j]=elevert[NV*i+j];
00083       }
00084       teleprop[k]=eleprop[i];
00085       k++;
00086     }
00087   }
00088   assert(k==n_ele-n_del);
00089 
00090   n_ele -= n_del;
00091 
00092   /* update counters */
00093   *pn_ele=n_ele;
00094 
00095   ierr = PetscFree(*pelevert);CHKERRQ(ierr);
00096   *pelevert=televert;
00097   ierr = PetscFree(*peleprop);CHKERRQ(ierr);
00098   *peleprop=teleprop;
00099 
00100   MagparFunctionLogReturn(0);
00101 }

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