1*17860365SKenneth E. Jansen #ifdef HAVE_PETSC
210167291SKenneth E. Jansen #include <stdio.h>
310167291SKenneth E. Jansen #include <stdlib.h>
410167291SKenneth E. Jansen #include <assert.h>
510167291SKenneth E. Jansen #include <math.h>
610167291SKenneth E. Jansen
710167291SKenneth E. Jansen #include <mpi.h>
810167291SKenneth E. Jansen #include <petsc.h>
910167291SKenneth E. Jansen
1010167291SKenneth E. Jansen #include "common_c.h"
1110167291SKenneth E. Jansen #include "FCMangle.h"
1210167291SKenneth E. Jansen
1310167291SKenneth E. Jansen #define fillsparsecpetscs FortranCInterface_GLOBAL_(fillsparsecpetscs, FILLSPARSECPETSCS)
1410167291SKenneth E. Jansen #define fillsparsecpetscc FortranCInterface_GLOBAL_(fillsparsecpetscc, FILLSPARSECPETSCC)
1510167291SKenneth E. Jansen
1610167291SKenneth E. Jansen #define COLMAJ2D(row,col,numrow) (row-1)+(col-1)*numrow
1710167291SKenneth E. Jansen #define COLMAJ3D(a,b,c,amax,bmax,cmax) (a-1)+amax*((b-1)+bmax*(c-1))
1810167291SKenneth E. Jansen #define ROWMAJ2D_ONE(row,col,numcol) (row-1)*numcol+(col-1)
1910167291SKenneth E. Jansen typedef long long int gcorp_t;
2010167291SKenneth E. Jansen
fillsparsecpetscs(gcorp_t * ieng,double * EGmass,Mat * lhsP)2110167291SKenneth E. Jansen void fillsparsecpetscs(gcorp_t* ieng, double* EGmass, Mat* lhsP)
2210167291SKenneth E. Jansen {
2310167291SKenneth E. Jansen int npro = propar.npro;
2410167291SKenneth E. Jansen int nshl = shpdat.nshl;
2510167291SKenneth E. Jansen double* mb = (double*) malloc(sizeof(double)*nshl*nshl); //block to insert
2610167291SKenneth E. Jansen int e,i,j,aa; //following along with fillsparse.f
2710167291SKenneth E. Jansen PetscInt* locat = (PetscInt*) malloc(sizeof(PetscInt)*nshl);
2810167291SKenneth E. Jansen for(e=0;e<npro;e++)
2910167291SKenneth E. Jansen {
3010167291SKenneth E. Jansen for(aa=0;aa<nshl;aa++) locat[aa]=ieng[e+npro*aa]-1;
3110167291SKenneth E. Jansen // for(aa=0;aa<nshl;aa++) assert(locat[aa]>=0);
3210167291SKenneth E. Jansen for (i=0; i<nshl; i++) { // fill up Ke with respective egmass
3310167291SKenneth E. Jansen for (j=0; j<nshl; j++) {
3410167291SKenneth E. Jansen mb[nshl*i + j] = EGmass[e + npro*(i + nshl*j)];
3510167291SKenneth E. Jansen }
3610167291SKenneth E. Jansen }
3710167291SKenneth E. Jansen //MatSetValuesBlocked(*lhsP, nshl , locat, nshl, locat, mb, ADD_VALUES);
3810167291SKenneth E. Jansen PetscInt petsc_nshl;
3910167291SKenneth E. Jansen petsc_nshl = (PetscInt) nshl;
4010167291SKenneth E. Jansen MatSetValues(*lhsP, petsc_nshl , locat, petsc_nshl, locat, mb, ADD_VALUES);
4110167291SKenneth E. Jansen }
4210167291SKenneth E. Jansen free(mb);
4310167291SKenneth E. Jansen free(locat);
4410167291SKenneth E. Jansen }
fillsparsecpetscc(gcorp_t * ieng,double * EGmass,Mat * lhsP)4510167291SKenneth E. Jansen void fillsparsecpetscc(gcorp_t* ieng, double* EGmass, Mat* lhsP)
4610167291SKenneth E. Jansen {
4710167291SKenneth E. Jansen int npro = propar.npro;
4810167291SKenneth E. Jansen int nshl = shpdat.nshl;
4910167291SKenneth E. Jansen int nedof = conpar.nedof;
5010167291SKenneth E. Jansen double* mb = (double*) malloc(sizeof(double)*nedof*nedof); //block to insert
5110167291SKenneth E. Jansen int e,i,j,aa; //following along with fillsparse.f
5210167291SKenneth E. Jansen //int* locat = (int*) malloc(sizeof(int)*nshl);
5310167291SKenneth E. Jansen PetscInt* locat = (PetscInt*) malloc(sizeof(PetscInt)*nshl);
5410167291SKenneth E. Jansen for(e=0;e<npro;e++)
5510167291SKenneth E. Jansen {
5610167291SKenneth E. Jansen for(aa=0;aa<nshl;aa++) locat[aa]=ieng[e+npro*aa]-1;
5710167291SKenneth E. Jansen // for(aa=0;aa<nshl;aa++) assert(locat[aa]>=0);
5810167291SKenneth E. Jansen for (i=0; i<nedof; i++) { /* fill up Ke with respective egmass */
5910167291SKenneth E. Jansen for (j=0; j<nedof; j++) {
6010167291SKenneth E. Jansen mb[nedof*i + j] = EGmass[e + npro*(i + nedof*j)];
6110167291SKenneth E. Jansen }
6210167291SKenneth E. Jansen }
6310167291SKenneth E. Jansen //MatSetValuesBlocked(*lhsP, nshl , locat, nshl, locat, mb, ADD_VALUES);
6410167291SKenneth E. Jansen PetscInt petsc_nshl;
6510167291SKenneth E. Jansen petsc_nshl = (PetscInt) nshl;
6610167291SKenneth E. Jansen MatSetValuesBlocked(*lhsP, petsc_nshl , locat, petsc_nshl, locat, mb, ADD_VALUES);
6710167291SKenneth E. Jansen }
6810167291SKenneth E. Jansen free(mb);
6910167291SKenneth E. Jansen free(locat);
7010167291SKenneth E. Jansen }
71*17860365SKenneth E. Jansen #endif
72*17860365SKenneth E. Jansen
73