xref: /phasta/phSolver/common/fillsparse.c (revision 1016729149754f57cd03fe576ba6fd0f1723ab31)
1*10167291SKenneth E. Jansen #include <stdio.h>
2*10167291SKenneth E. Jansen #include <stdlib.h>
3*10167291SKenneth E. Jansen #include <assert.h>
4*10167291SKenneth E. Jansen #include <math.h>
5*10167291SKenneth E. Jansen 
6*10167291SKenneth E. Jansen #include <mpi.h>
7*10167291SKenneth E. Jansen #include <petsc.h>
8*10167291SKenneth E. Jansen 
9*10167291SKenneth E. Jansen #include "common_c.h"
10*10167291SKenneth E. Jansen #include "FCMangle.h"
11*10167291SKenneth E. Jansen 
12*10167291SKenneth E. Jansen #define fillsparsecpetscs FortranCInterface_GLOBAL_(fillsparsecpetscs, FILLSPARSECPETSCS)
13*10167291SKenneth E. Jansen #define fillsparsecpetscc FortranCInterface_GLOBAL_(fillsparsecpetscc, FILLSPARSECPETSCC)
14*10167291SKenneth E. Jansen 
15*10167291SKenneth E. Jansen #define COLMAJ2D(row,col,numrow) (row-1)+(col-1)*numrow
16*10167291SKenneth E. Jansen #define COLMAJ3D(a,b,c,amax,bmax,cmax) (a-1)+amax*((b-1)+bmax*(c-1))
17*10167291SKenneth E. Jansen #define ROWMAJ2D_ONE(row,col,numcol) (row-1)*numcol+(col-1)
18*10167291SKenneth E. Jansen typedef long long int gcorp_t;
19*10167291SKenneth E. Jansen 
20*10167291SKenneth E. Jansen void fillsparsecpetscs(gcorp_t* ieng, double* EGmass, Mat* lhsP)
21*10167291SKenneth E. Jansen {
22*10167291SKenneth E. Jansen         int npro = propar.npro;
23*10167291SKenneth E. Jansen         int nshl = shpdat.nshl;
24*10167291SKenneth E. Jansen         double* mb = (double*) malloc(sizeof(double)*nshl*nshl); //block to insert
25*10167291SKenneth E. Jansen         int e,i,j,aa; //following along with fillsparse.f
26*10167291SKenneth E. Jansen         PetscInt* locat = (PetscInt*) malloc(sizeof(PetscInt)*nshl);
27*10167291SKenneth E. Jansen         for(e=0;e<npro;e++)
28*10167291SKenneth E. Jansen         {
29*10167291SKenneth E. Jansen          for(aa=0;aa<nshl;aa++) locat[aa]=ieng[e+npro*aa]-1;
30*10167291SKenneth E. Jansen //         for(aa=0;aa<nshl;aa++) assert(locat[aa]>=0);
31*10167291SKenneth E. Jansen          for (i=0; i<nshl; i++)  {   // fill up Ke with respective egmass
32*10167291SKenneth E. Jansen            for (j=0; j<nshl; j++)  {
33*10167291SKenneth E. Jansen             mb[nshl*i + j] = EGmass[e + npro*(i + nshl*j)];
34*10167291SKenneth E. Jansen            }
35*10167291SKenneth E. Jansen          }
36*10167291SKenneth E. Jansen          //MatSetValuesBlocked(*lhsP, nshl , locat, nshl, locat, mb, ADD_VALUES);
37*10167291SKenneth E. Jansen          PetscInt petsc_nshl;
38*10167291SKenneth E. Jansen          petsc_nshl = (PetscInt) nshl;
39*10167291SKenneth E. Jansen          MatSetValues(*lhsP, petsc_nshl , locat, petsc_nshl, locat, mb, ADD_VALUES);
40*10167291SKenneth E. Jansen         }
41*10167291SKenneth E. Jansen         free(mb);
42*10167291SKenneth E. Jansen 	free(locat);
43*10167291SKenneth E. Jansen }
44*10167291SKenneth E. Jansen void fillsparsecpetscc(gcorp_t* ieng, double* EGmass, Mat* lhsP)
45*10167291SKenneth E. Jansen {
46*10167291SKenneth E. Jansen         int npro = propar.npro;
47*10167291SKenneth E. Jansen         int nshl = shpdat.nshl;
48*10167291SKenneth E. Jansen         int nedof = conpar.nedof;
49*10167291SKenneth E. Jansen         double* mb = (double*) malloc(sizeof(double)*nedof*nedof); //block to insert
50*10167291SKenneth E. Jansen         int e,i,j,aa; //following along with fillsparse.f
51*10167291SKenneth E. Jansen         //int* locat = (int*) malloc(sizeof(int)*nshl);
52*10167291SKenneth E. Jansen         PetscInt* locat = (PetscInt*) malloc(sizeof(PetscInt)*nshl);
53*10167291SKenneth E. Jansen         for(e=0;e<npro;e++)
54*10167291SKenneth E. Jansen         {
55*10167291SKenneth E. Jansen          for(aa=0;aa<nshl;aa++) locat[aa]=ieng[e+npro*aa]-1;
56*10167291SKenneth E. Jansen //         for(aa=0;aa<nshl;aa++) assert(locat[aa]>=0);
57*10167291SKenneth E. Jansen          for (i=0; i<nedof; i++)  {   /* fill up Ke with respective egmass */
58*10167291SKenneth E. Jansen            for (j=0; j<nedof; j++)  {
59*10167291SKenneth E. Jansen             mb[nedof*i + j] = EGmass[e + npro*(i + nedof*j)];
60*10167291SKenneth E. Jansen            }
61*10167291SKenneth E. Jansen          }
62*10167291SKenneth E. Jansen          //MatSetValuesBlocked(*lhsP, nshl , locat, nshl, locat, mb, ADD_VALUES);
63*10167291SKenneth E. Jansen          PetscInt petsc_nshl;
64*10167291SKenneth E. Jansen          petsc_nshl = (PetscInt) nshl;
65*10167291SKenneth E. Jansen          MatSetValuesBlocked(*lhsP, petsc_nshl , locat, petsc_nshl, locat, mb, ADD_VALUES);
66*10167291SKenneth E. Jansen         }
67*10167291SKenneth E. Jansen         free(mb);
68*10167291SKenneth E. Jansen 	free(locat);
69*10167291SKenneth E. Jansen }
70