xref: /phasta/phSolver/incompressible/usr.c (revision 18924a0afcafce200c2f5b454f3c68488e192bf0)
1 /*===========================================================================
2  *
3  * "usr.c":  user's function
4  *
5  *===========================================================================
6  */
7 #include <stdio.h>
8 #include <stdlib.h>
9 #include "les.h"
10 #include "usr.h"
11 #include "common_c.h"
12 #include "phastaIO.h"
13 #include "phIO.h"
14 #include "syncio.h"
15 #include "posixio.h"
16 #include "streamio.h"
17 #include "new_interface.h"
18 #include <FCMangle.h>
19 
20 extern char phasta_iotype[80];
21 
22 /*===========================================================================
23  *
24  * "usrNew":  Put all the values in usrHd
25  *
26  * From FORTRAN
27  *
28  *	integer		usr(100)
29  *	dimension	aperm(numnp,nperm)
30  *	...
31  *	call usrnew( usr, aperm, ..., numnp, ...)
32  *
33  *
34  *===========================================================================
35  */
36 #include "mpi.h"
37 static int lmNum = 0;
38 static LesHd lesArray[8];
39 void   usrNew(	UsrHd	  usrHd,
40                         int*      eqnType,
41                         double*	  aperm,
42                         double*	  atemp,
43                         double*   resf,
44                         double*   solinc,
45                         double*   flowDiag,
46                         double*   sclrDiag,
47                         double*   lesP,
48                         double*   lesQ,
49                         int*      iBC,
50                         double*   BC,
51                         int*      iper,
52                         int*      ilwork,
53                         int*      numpe,
54                         int*      nNodes,
55                         int*      nenl,
56                         int*	  nPermDims,
57                         int*	  nTmpDims,
58                         int*	  rowp,
59                         int*	  colm,
60                         double*   lhsK,
61                         double*   lhsP,
62                         double*   lhsS,
63                         int*      nnz_tot,
64                         double*   CGsol
65     )
66 {
67     char*	funcName = "usrNew" ;	/* function name		*/
68 
69 /*---------------------------------------------------------------------------
70  * Stick the parameters
71  *---------------------------------------------------------------------------
72  */
73     usrHd->eqnType      = *eqnType ;
74     usrHd->aperm	= aperm ;
75     usrHd->atemp	= atemp ;
76     usrHd->resf         = resf ;
77     usrHd->solinc       = solinc ;
78     usrHd->flowDiag     = flowDiag ;
79     usrHd->sclrDiag     = sclrDiag ;
80     usrHd->lesP         = lesP ;
81     usrHd->lesQ         = lesQ ;
82     usrHd->iBC          = iBC  ;
83     usrHd->BC           = BC   ;
84     usrHd->iper         = iper ;
85     usrHd->ilwork       = ilwork ;
86     usrHd->numpe        = *numpe ;
87     usrHd->nNodes	= *nNodes ;
88     usrHd->nenl         = *nenl ;
89     usrHd->nPermDims	= *nPermDims ;
90     usrHd->nTmpDims	= *nTmpDims ;
91     usrHd->rowp	        = rowp ;
92     usrHd->colm	        = colm ;
93     usrHd->lhsK	        = lhsK ;
94     usrHd->lhsP	        = lhsP ;
95     usrHd->lhsS         = lhsS ;
96     usrHd->nnz_tot      = nnz_tot ;
97     usrHd->CGsol        = CGsol;
98 } /* end of usrNew() */
99 
100 /*===========================================================================
101  *
102  * "usrPointer":  Get the pointer
103  *
104  *===========================================================================
105  */
106 Real*
107 usrPointer(	UsrHd	usrHd,
108             Integer	id,
109             Integer	offset,
110             Integer	nDims )
111 {
112     char*	funcName = "usrPointer";/* function name		*/
113     Real*	pnt ;			/* pointer			*/
114 
115 /*---------------------------------------------------------------------------
116  * Get the head of the memory
117  *---------------------------------------------------------------------------
118  */
119     if ( id == LES_RES_PNT ) {
120 
121         pnt	= usrHd->resf ;
122         id	= 0 ;
123 
124     } else if ( id == LES_SOL_PNT ) {
125 
126         pnt	= usrHd->solinc ;
127         id	= 0 ;
128 
129     } else if ( id < 0 ) {
130 
131         pnt	= usrHd->aperm ;
132         id	= id + usrHd->nPermDims ;
133 
134     } else {
135 
136         pnt	= usrHd->atemp ;
137         id	= id ;
138 
139     }
140 /*---------------------------------------------------------------------------
141  * Get the offset
142  *---------------------------------------------------------------------------
143  */
144     pnt		= pnt + (id + offset) * usrHd->nNodes ;
145 
146 /*---------------------------------------------------------------------------
147  * Return the pointer
148  *---------------------------------------------------------------------------
149  */
150     return( pnt ) ;
151 
152 } /* end of usrPointer() */
153 
154 #define myflesnew FortranCInterface_GLOBAL_(myflesnew,MYFLESNEW)
155 #define myflessolve FortranCInterface_GLOBAL_(myflessolve,MYFLESSOLVE)
156 #define savelesrestart FortranCInterface_GLOBAL_(savelesrestart,SAVELESRESTART)
157 #define readlesrestart FortranCInterface_GLOBAL_(readlesrestart,READLESRESTART)
158 #define solverlicenseserver FortranCInterface_GLOBAL_(solverlicenseserver,SOLVERLICENSESERVER)
159 
160 
161 
162 #ifdef intel
163         lesArray[ *lesId ] = lesNew( fileName, *lmport, &lmNum, *eqnType,
164                                      *nDofs, *minIters, *maxIters, *nKvecs,
165                                      *prjFlag, *nPrjs, *presPrjFlag, *nPresPrjs,presPrec,
166                                      *tol, *presTol, *verbose, stats, nPermDims,
167                                      nTmpDims );
168     return ;}
169 /* the following is a fake function that was required when we moved to
170    a C++ main on in the MS Visual Studio environment.  It fails to
171    link because it is looking for this function
172 */
173 void  _CrtDbgReport() {
174     return ;}
175 
176 double __vcos_(double fg) { fflush(stdout); printf(" vcos got called \n"); fflush(stdout);}
177 double __vlog_(double fg)  { fflush(stdout); printf(" vlog got called \n"); fflush(stdout);}
178 
179 
180 #endif /* we are in unix land... whew.  secretly we have equivalenced fileName and  */
181 
182 /* #ifdef LINUX*/
183 /* void flush_(int* junk ){ return; }*/
184 /* #endif*/
185 void    myflesnew(	     Integer*	lesId,
186                          Integer*	lmport,
187                          Integer*	eqnType,
188                          Integer*	nDofs,
189                          Integer*	minIters,
190                          Integer*	maxIters,
191                          Integer*	nKvecs,
192                          Integer*	prjFlag,
193                          Integer*	nPrjs,
194                          Integer*	presPrjFlag,
195                          Integer*	nPresPrjs,
196                          Real*	    tol,
197                          Real*     	presTol,
198                          Integer*	verbose,
199                          Real*     	stats,
200                          Integer*	nPermDims,
201                          Integer*	nTmpDims,
202                          char*      lmhost          ) {
203     int procId;
204 #ifdef AMG
205     int presPrec=1;
206 #else
207     int presPrec=0;
208 #endif
209     MPI_Comm_rank( MPI_COMM_WORLD, &procId ) ;
210     if(lmNum==0){
211         if(procId==0){
212             lesArray[ *lesId ] = lesNew( lmhost, *lmport, &lmNum, *eqnType,
213                                          *nDofs, *minIters, *maxIters, *nKvecs,
214                                          *prjFlag, *nPrjs, *presPrjFlag, *nPresPrjs,presPrec,
215                                          *tol, *presTol, *verbose, stats, nPermDims,
216                                          nTmpDims );
217             MPI_Bcast( &lmNum, 1, MPI_INT, 0, MPI_COMM_WORLD ) ;
218         } else {
219             MPI_Bcast( &lmNum, 1, MPI_INT, 0, MPI_COMM_WORLD ) ;
220             lesArray[ *lesId ] = lesNew( lmhost, *lmport, &lmNum, *eqnType,
221                                          *nDofs, *minIters, *maxIters, *nKvecs,
222                                          *prjFlag, *nPrjs, *presPrjFlag, *nPresPrjs,presPrec,
223                                          *tol, *presTol, *verbose, stats, nPermDims,
224                                          nTmpDims );
225         }
226     } else {
227         lesArray[ *lesId ] = lesNew( lmhost, *lmport, &lmNum, *eqnType,
228                                      *nDofs, *minIters, *maxIters, *nKvecs,
229                                      *prjFlag, *nPrjs, *presPrjFlag, *nPresPrjs,presPrec,
230                                      *tol, *presTol, *verbose, stats, nPermDims,
231                                      nTmpDims );
232     }
233     return ;
234 }
235 
236 
237 void
238 savelesrestart( Integer* lesId,
239                  Real*    aperm,
240                  Integer* nshg,
241                  Integer* myrank,
242                  Integer* lstep,
243                  Integer* nPermDims ) {
244 
245     int nPrjs, PrjSrcId;
246     int nPresPrjs, PresPrjSrcId;
247     char filename[255];
248     int iarray[3];
249     int size, nitems;
250     double* projVec;
251     int i, j, count;
252 
253     nPrjs = (Integer) lesGetPar( lesArray[ *lesId ], LES_ACT_PRJS );
254     PrjSrcId = (Integer) lesGetPar( lesArray[ *lesId ], LES_PRJ_VEC_ID );
255 
256     if ( PrjSrcId < 0 ) PrjSrcId += *nPermDims;
257 
258     projVec = (double*)malloc( nPrjs * ( *nshg ) * sizeof( double ) );
259 
260     count = 0;
261     for( i = PrjSrcId; i < PrjSrcId+nPrjs; i ++ ) {
262         for( j = 0 ; j < *nshg; j++ ) {
263             projVec[ count++ ] = aperm[ (*nshg) * i + j ];
264         }
265     }
266 
267     iarray[ 0 ] = *nshg;
268     iarray[ 1 ] = nPrjs;
269     nitems = 2;
270     size = (*nshg)*nPrjs;
271 
272     int name_length;
273     name_length = 18;
274     Write_Field(myrank,"a","projection vectors",&name_length, (void *)projVec,"d", nshg, &nPrjs, lstep);
275 
276     free(projVec);
277 
278     nPresPrjs = (Integer) lesGetPar( lesArray[ *lesId ], LES_ACT_PRES_PRJS );
279     PresPrjSrcId =(Integer)lesGetPar( lesArray[ *lesId ], LES_PRES_PRJ_VEC_ID );
280     if ( PresPrjSrcId < 0 ) PresPrjSrcId += *nPermDims;
281 
282     projVec = (double*)malloc( nPresPrjs * ( *nshg ) * sizeof( double ) );
283 
284     count = 0;
285     for( i = PresPrjSrcId; i < (PresPrjSrcId + nPresPrjs) ; i ++ ) {
286         for( j = 0 ; j < *nshg; j++ ) {
287             projVec[ count++ ] = aperm[ (*nshg) * i + j ];
288         }
289     }
290 
291     iarray[ 0 ] = *nshg;
292     iarray[ 1 ] = nPresPrjs;
293     nitems = 2;
294     size = (*nshg)*nPresPrjs;
295 
296     name_length = 27;
297     Write_Field(myrank,"a","pressure projection vectors",&name_length, projVec,"d", nshg, &nPresPrjs, lstep);
298 
299     free( projVec);
300 }
301 
302 void
303 readlesrestart( Integer* lesId,
304                  Real*    aperm,
305                  Integer* nshg,
306                  Integer* myrank,
307                  Integer* lstep ,
308                  Integer* nPermDims ) {
309 
310     int nPrjs, PrjSrcId;
311     int nPresPrjs, PresPrjSrcId;
312     char filename[255];
313     phio_fp fileHandle = NULL;
314     int iarray[3]={-1,-1,-1};
315     int size, nitems;
316     int itwo=2;
317     int lnshg;
318     double* projVec;
319     int i,j,count;
320     int nfiles;
321     int nfields;
322     int numParts;
323     int nprocs;
324     int nppf;
325 
326     nfiles = outpar.nsynciofiles;
327     numParts = workfc.numpe;
328     nprocs = workfc.numpe;
329     // Calculate number of parts each proc deal with and where it start and end ...
330     int nppp = numParts/nprocs;        // nppp : Number of parts per proc ...
331     int startpart = *myrank * nppp +1;    // Part id from which I (myrank) start ...
332     int endpart = startpart + nppp - 1;  // Part id to which I (myrank) end ...
333 
334     if( nfiles == -1 )
335       streamio_setup_read(&fileHandle, streamio_get_gr());
336     else if( nfiles == 0 )
337       posixio_setup(&fileHandle, 'r');
338     else if( nfiles > 0 )
339       syncio_setup_read(nfiles, &fileHandle);
340     phio_constructName(fileHandle,"restart",filename);
341     phio_appendInt(filename, *lstep);
342     phio_openfile(filename, fileHandle);
343 
344     if ( !fileHandle ) return; // See phastaIO.cc for error fileHandle
345     phio_readheader(fileHandle, "projection vectors", (void*)iarray,
346                 &itwo, "integer", phasta_iotype);
347 
348     if ( iarray[0] != *nshg ) {
349         phio_closefile(fileHandle);
350         if(workfc.myrank==workfc.master)
351           printf("projection vectors are being initialized to zero (SAFE)\n");
352         return;
353     }
354 
355     lnshg = iarray[ 0 ] ;
356     nPrjs = iarray[ 1 ] ;
357 
358     size = (*nshg)*nPrjs;
359     projVec = (double*)malloc( size * sizeof( double ));
360 
361     phio_readdatablock(fileHandle, "projection vectors", (void*)projVec,
362                     &size, "double", phasta_iotype );
363 
364     lesSetPar( lesArray[ *lesId ], LES_ACT_PRJS, (Real) nPrjs );
365     PrjSrcId = (Integer) lesGetPar( lesArray[ *lesId ], LES_PRJ_VEC_ID );
366     if ( PrjSrcId < 0 ) PrjSrcId += *nPermDims;
367 
368     count = 0;
369     for( i = PrjSrcId; i < PrjSrcId+nPrjs; i ++ ) {
370         for( j = 0 ; j < *nshg; j++ ) {
371             aperm[ (*nshg) * i + j ] = projVec[ count++ ] ;
372         }
373     }
374 
375     free( projVec );
376 
377     iarray[0] = -1; iarray[1] = -1; iarray[2] = -1;
378 
379     phio_readheader(fileHandle, "pressure projection vectors", (void*)iarray,
380                  &itwo, "integer", phasta_iotype );
381 
382     lnshg = iarray[ 0 ] ;
383     nPresPrjs = iarray[ 1 ] ;
384 
385     if ( lnshg != *nshg )  {
386         phio_closefile(fileHandle);
387         if(workfc.myrank==workfc.master)
388           printf("pressure projection vectors are being initialized to zero (SAFE)\n");
389         return;
390     }
391 
392     size = (*nshg)*nPresPrjs;
393     projVec = (double*)malloc( size * sizeof( double ));
394 
395     phio_readdatablock(fileHandle, "pressure projection vectors", (void*)projVec,
396                     &size, "double", phasta_iotype );
397 
398     lesSetPar( lesArray[ *lesId ], LES_ACT_PRES_PRJS, (Real) nPresPrjs );
399     PresPrjSrcId=(Integer)lesGetPar( lesArray[ *lesId ], LES_PRES_PRJ_VEC_ID );
400     if ( PresPrjSrcId < 0 ) PresPrjSrcId += *nPermDims;
401 
402     count = 0;
403     for( i = PresPrjSrcId; i < PresPrjSrcId+nPresPrjs; i ++ ) {
404         for( j = 0 ; j < *nshg; j++ ) {
405             aperm[ (*nshg) * i + j ] = projVec[ count++ ] ;
406         }
407     }
408 
409     free( projVec );
410 
411     phio_closefile(fileHandle);
412 }
413 
414 void  myflessolve( Integer* lesId,
415                     UsrHd    usrHd){
416     lesSolve( lesArray[ *lesId ], usrHd );
417 }
418 
419 
420 int solverlicenseserver(char key[]){
421 #ifdef intel
422     strcpy(key,"C:\\cygwin\\license.dat");
423 #else
424     char* env_server_name;
425     env_server_name = getenv("LES_LICENSE_SERVER");
426     if(env_server_name) strcpy(key, env_server_name);
427     else {
428         if(workfc.myrank==workfc.master) {
429           fprintf(stderr,"environment variable LES_LICENSE_SERVER not defined \n");
430           fprintf(stderr,"using wesley as default \n");
431         }
432         strcpy(key, "acusim.license.scorec.rpi.edu");
433     }
434 #endif
435     return 1;
436 }
437