Ice Sheet System Model  4.18
Code documentation
All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros
Functions
MpiDenseMumpsSolve.cpp File Reference

: solve dense matrix system with MUMPS More...

#include "../../shared/Numerics/types.h"
#include "../../shared/MemOps/MemOps.h"
#include "../../shared/Exceptions/exceptions.h"
#include "../../shared/io/Comm/IssmComm.h"
#include "../../classes/Params/GenericParam.h"
#include "../../classes/Params/Parameters.h"
#include "../mpi/issmmpi.h"
#include "../adolc/adolcincludes.h"
#include "./mumpsincludes.h"
#include <dmumps_c.h>

Go to the source code of this file.

Functions

void MumpsInit (DMUMPS_STRUC_C &theMumpsStruc)
 
void MumpsSettings (DMUMPS_STRUC_C &theMumpsStruc)
 
void MumpsAnalyze (DMUMPS_STRUC_C &theMumpsStruc)
 
void MumpsFactorize (DMUMPS_STRUC_C &theMumpsStruc)
 
void MumpsBacksubstitute (DMUMPS_STRUC_C &theMumpsStruc)
 
void MumpsFinalize (DMUMPS_STRUC_C &theMumpsStruc)
 
void MumpsSolve (int n, int nnz, int local_nnz, int *irn_loc, int *jcn_loc, IssmPDouble *a_loc, IssmPDouble *rhs, Parameters *parameters=0)
 
void MpiDenseMumpsSolve (IssmDouble *uf, int uf_M, int uf_m, IssmDouble *Kff, int Kff_M, int Kff_N, int Kff_m, IssmDouble *pf, int pf_M, int pf_m, Parameters *parameters)
 

Detailed Description

: solve dense matrix system with MUMPS

Definition in file MpiDenseMumpsSolve.cpp.

Function Documentation

◆ MumpsInit()

void MumpsInit ( DMUMPS_STRUC_C &  theMumpsStruc)

Definition at line 25 of file MpiDenseMumpsSolve.cpp.

25  {
26  theMumpsStruc.par = 1;
27  theMumpsStruc.sym = 0;
28  theMumpsStruc.comm_fortran = MPI_Comm_c2f(IssmComm::GetComm());
29  theMumpsStruc.job = -1;
30  dmumps_c(&theMumpsStruc);
31 }

◆ MumpsSettings()

void MumpsSettings ( DMUMPS_STRUC_C &  theMumpsStruc)

Definition at line 34 of file MpiDenseMumpsSolve.cpp.

34  {
35  /*Control statements:{{{ */
36  theMumpsStruc.icntl[1-1] = 6; //error verbose: default 6, 0 or negative -> suppressed
37  theMumpsStruc.icntl[2-1] = 0; //std verbose: default 1, 0 or negative -> suppressed
38  theMumpsStruc.icntl[3-1] = 0; //global information verbose: default 6, 0 or negative -> suppressed
39  theMumpsStruc.icntl[4-1] = 0; //verbose everything: default is 4
40  theMumpsStruc.icntl[5-1] = 0;
41  theMumpsStruc.icntl[18-1] = 3;
42 
43  theMumpsStruc.icntl[20-1] = 0;
44  theMumpsStruc.icntl[21-1] = 0;
45  theMumpsStruc.icntl[30-1] = 0;
46  /*}}}*/
47 }

◆ MumpsAnalyze()

void MumpsAnalyze ( DMUMPS_STRUC_C &  theMumpsStruc)

Definition at line 50 of file MpiDenseMumpsSolve.cpp.

50  {
51  theMumpsStruc.job = 1;
52  dmumps_c(&theMumpsStruc);
53 }

◆ MumpsFactorize()

void MumpsFactorize ( DMUMPS_STRUC_C &  theMumpsStruc)

Definition at line 56 of file MpiDenseMumpsSolve.cpp.

56  {
57  theMumpsStruc.job = 2;
58  dmumps_c(&theMumpsStruc);
59 }

◆ MumpsBacksubstitute()

void MumpsBacksubstitute ( DMUMPS_STRUC_C &  theMumpsStruc)

Definition at line 62 of file MpiDenseMumpsSolve.cpp.

62  {
63  theMumpsStruc.job = 3;
64  dmumps_c(&theMumpsStruc);
65 }

◆ MumpsFinalize()

void MumpsFinalize ( DMUMPS_STRUC_C &  theMumpsStruc)

Definition at line 68 of file MpiDenseMumpsSolve.cpp.

68  {
69  theMumpsStruc.job = -2;
70  dmumps_c(&theMumpsStruc);
71 }

◆ MumpsSolve()

void MumpsSolve ( int  n,
int  nnz,
int  local_nnz,
int *  irn_loc,
int *  jcn_loc,
IssmPDouble a_loc,
IssmPDouble rhs,
Parameters parameters = 0 
)

Definition at line 73 of file MpiDenseMumpsSolve.cpp.

73  {
74  /*Initialize mumps*/
75  DMUMPS_STRUC_C theMumpsStruc;
76  MumpsInit(theMumpsStruc);
77  MumpsSettings(theMumpsStruc);
78 
79  /*now setup the rest of theMumpsStruc */
80  theMumpsStruc.n=n;
81  theMumpsStruc.nz=nnz;
82  theMumpsStruc.nz_loc=local_nnz;
83  theMumpsStruc.irn_loc=irn_loc;
84  theMumpsStruc.jcn_loc=jcn_loc;
85  theMumpsStruc.a_loc=a_loc;
86  theMumpsStruc.rhs=rhs;
87  theMumpsStruc.nrhs=1;
88  theMumpsStruc.lrhs=1;
89 
90  /*Solve system*/
91  MumpsAnalyze(theMumpsStruc);
92  MumpsFactorize(theMumpsStruc);
93  MumpsBacksubstitute(theMumpsStruc);
94  MumpsFinalize(theMumpsStruc);
95 }

◆ MpiDenseMumpsSolve()

void MpiDenseMumpsSolve ( IssmDouble uf,
int  uf_M,
int  uf_m,
IssmDouble Kff,
int  Kff_M,
int  Kff_N,
int  Kff_m,
IssmDouble pf,
int  pf_M,
int  pf_m,
Parameters parameters 
)

Definition at line 109 of file MpiDenseMumpsSolve.cpp.

109  { /*{{{*/
110 
111  /*Variables*/
112  ISSM_MPI_Comm comm;
113  int num_procs;
114  int i,j;
115  int nnz,local_nnz;
116  int *irn_loc = NULL;
117  int *jcn_loc = NULL;
118  IssmDouble *a_loc = NULL;
119  int count;
120  int lower_row;
121  int upper_row;
122  IssmDouble *rhs = NULL;
123  int *recvcounts = NULL;
124  int *displs = NULL;
125 
126  /*Communicator info */
127  num_procs = IssmComm::GetSize();
128  comm = IssmComm::GetComm();
129 
130  /*First, some checks*/
131  if (Kff_M!=Kff_N)_error_("stiffness matrix Kff should be square");
132  if (uf_M!=Kff_M || uf_M!=pf_M)_error_("solution vector should be the same size as stiffness matrix Kff and load vector pf");
133  if (uf_m!=Kff_m || uf_m!=pf_m)_error_("solution vector should be locally the same size as stiffness matrix Kff and load vector pf");
134 
135  /*Initialize matrix */
136  /*figure out number of non-zero entries: */
137  local_nnz=0;
138  for(i=0;i<Kff_m;i++){
139  for(j=0;j<Kff_N;j++){
140  if (Kff[i*Kff_N+j]!=0)local_nnz++;
141  }
142  }
143 
144  ISSM_MPI_Reduce(&local_nnz,&nnz,1,ISSM_MPI_INT,ISSM_MPI_SUM,0,comm);
145  ISSM_MPI_Bcast(&nnz,1,ISSM_MPI_INT,0,comm);
146 
147  /*Allocate: */
148  if(local_nnz){
149  irn_loc=xNew<int>(local_nnz);
150  jcn_loc=xNew<int>(local_nnz);
151  a_loc=xNew<IssmDouble>(local_nnz);
152  }
153 
154  /*Populate the triplets: */
155  GetOwnershipBoundariesFromRange(&lower_row,&upper_row,Kff_m,comm);
156  count=0;
157  for(i=0;i<Kff_m;i++){
158  for(j=0;j<Kff_N;j++){
159  if (Kff[i*Kff_N+j]!=0){
160  irn_loc[count]=lower_row+i+1; //fortran indexing
161  jcn_loc[count]=j+1; //fortran indexing
162  a_loc[count]=Kff[i*Kff_N+j];
163  count++;
164  }
165  }
166  }
167  /*Deal with right hand side. We need to ISSM_MPI_Gather it onto cpu 0: */
168  rhs=xNew<IssmDouble>(pf_M);
169 
170  recvcounts=xNew<int>(num_procs);
171  displs=xNew<int>(num_procs);
172 
173  /*recvcounts:*/
174  ISSM_MPI_Allgather(&pf_m,1,ISSM_MPI_INT,recvcounts,1,ISSM_MPI_INT,comm);
175 
176  /*displs: */
177  ISSM_MPI_Allgather(&lower_row,1,ISSM_MPI_INT,displs,1,ISSM_MPI_INT,comm);
178 
179  /*Gather:*/
180  ISSM_MPI_Gatherv(pf, pf_m, ISSM_MPI_DOUBLE, rhs, recvcounts, displs, ISSM_MPI_DOUBLE,0,comm);
181 
182  MumpsSolve(Kff_M,nnz,local_nnz,irn_loc,jcn_loc,a_loc,rhs,parameters);
183 
184  /*Now scatter from cpu 0 to all other cpus*/
185  ISSM_MPI_Scatterv( rhs, recvcounts, displs, ISSM_MPI_DOUBLE, uf, uf_m, ISSM_MPI_DOUBLE, 0, comm);
186 
187  /*Cleanup*/
188  xDelete<int>(irn_loc);
189  xDelete<int>(jcn_loc);
190  xDelete<IssmDouble>(a_loc);
191  xDelete<IssmDouble>(rhs);
192  xDelete<int>(recvcounts);
193  xDelete<int>(displs);
194 } /*}}}*/
IssmDouble
double IssmDouble
Definition: types.h:37
MumpsAnalyze
void MumpsAnalyze(DMUMPS_STRUC_C &theMumpsStruc)
Definition: MpiDenseMumpsSolve.cpp:50
ISSM_MPI_Gatherv
int ISSM_MPI_Gatherv(void *sendbuf, int sendcnt, ISSM_MPI_Datatype sendtype, void *recvbuf, int *recvcnts, int *displs, ISSM_MPI_Datatype recvtype, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:279
ISSM_MPI_SUM
#define ISSM_MPI_SUM
Definition: issmmpi.h:134
IssmComm::GetComm
static ISSM_MPI_Comm GetComm(void)
Definition: IssmComm.cpp:30
GetOwnershipBoundariesFromRange
void GetOwnershipBoundariesFromRange(int *plower_row, int *pupper_row, int range, ISSM_MPI_Comm comm)
Definition: GetOwnershipBoundariesFromRange.cpp:16
ISSM_MPI_DOUBLE
#define ISSM_MPI_DOUBLE
Definition: issmmpi.h:125
ISSM_MPI_INT
#define ISSM_MPI_INT
Definition: issmmpi.h:127
IssmComm::GetSize
static int GetSize(void)
Definition: IssmComm.cpp:46
MumpsInit
void MumpsInit(DMUMPS_STRUC_C &theMumpsStruc)
Definition: MpiDenseMumpsSolve.cpp:25
ISSM_MPI_Bcast
int ISSM_MPI_Bcast(void *buffer, int count, ISSM_MPI_Datatype datatype, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:162
ISSM_MPI_Allgather
int ISSM_MPI_Allgather(void *sendbuf, int sendcount, ISSM_MPI_Datatype sendtype, void *recvbuf, int recvcount, ISSM_MPI_Datatype recvtype, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:45
_error_
#define _error_(StreamArgs)
Definition: exceptions.h:49
ISSM_MPI_Comm
int ISSM_MPI_Comm
Definition: issmmpi.h:118
MumpsBacksubstitute
void MumpsBacksubstitute(DMUMPS_STRUC_C &theMumpsStruc)
Definition: MpiDenseMumpsSolve.cpp:62
MumpsFinalize
void MumpsFinalize(DMUMPS_STRUC_C &theMumpsStruc)
Definition: MpiDenseMumpsSolve.cpp:68
MumpsSettings
void MumpsSettings(DMUMPS_STRUC_C &theMumpsStruc)
Definition: MpiDenseMumpsSolve.cpp:34
ISSM_MPI_Reduce
int ISSM_MPI_Reduce(void *sendbuf, void *recvbuf, int count, ISSM_MPI_Datatype datatype, ISSM_MPI_Op op, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:373
MumpsFactorize
void MumpsFactorize(DMUMPS_STRUC_C &theMumpsStruc)
Definition: MpiDenseMumpsSolve.cpp:56
MumpsSolve
void MumpsSolve(int n, int nnz, int local_nnz, int *irn_loc, int *jcn_loc, IssmPDouble *a_loc, IssmPDouble *rhs, Parameters *parameters=0)
Definition: MpiDenseMumpsSolve.cpp:73
ISSM_MPI_Scatterv
int ISSM_MPI_Scatterv(void *sendbuf, int *sendcnts, int *displs, ISSM_MPI_Datatype sendtype, void *recvbuf, int recvcnt, ISSM_MPI_Datatype recvtype, int root, ISSM_MPI_Comm comm)
Definition: issmmpi.cpp:444